Actual source code: snesngmres.h
1: #pragma once
3: #include <petsc/private/snesimpl.h>
5: /* Data structure for the Nonlinear GMRES method. */
6: typedef struct {
7: /* Solver parameters and counters */
8: PetscInt msize; /* maximum size of krylov space */
9: PetscInt restart_it; /* number of iterations the restart conditions persist before restart */
10: PetscViewer monitor; /* debugging output for NGMRES */
11: PetscInt restart_periodic; /* number of iterations to restart after */
13: SNESNGMRESRestartType restart_type;
14: SNESNGMRESSelectType select_type;
16: /* History and subspace data */
17: Vec *Fdot; /* residual history -- length msize */
18: Vec *Xdot; /* solution history -- length msize */
19: PetscReal *fnorms; /* the residual norm history */
20: PetscReal *xnorms; /* the solution norm history */
22: /* General minimization problem context */
23: PetscScalar *h; /* the constraint matrix */
24: PetscScalar *beta; /* rhs for the minimization problem */
25: PetscScalar *xi; /* the dot-product of the current and previous res. */
27: /* Line searches */
28: SNESLineSearch additive_linesearch; /* Line search for the additive variant */
30: /* Selection constants */
31: PetscBool candidate; /* use candidate storage approach */
32: PetscBool approxfunc; /* approximate the function rather than recomputing it */
33: PetscBool singlereduction; /* use a single reduction (with more local work) for tolerance selection */
34: PetscReal gammaA; /* Criterion A residual tolerance */
35: PetscReal epsilonB; /* Criterion B difference tolerance */
36: PetscReal deltaB; /* Criterion B residual tolerance */
37: PetscReal gammaC; /* Restart residual tolerance */
38: PetscBool restart_fm_rise; /* Restart on F_M residual increase */
40: PetscReal andersonBeta; /* Relaxation parameter for Anderson Mixing */
42: /* Least squares minimization solve context */
43: PetscScalar *q; /* the matrix formed as q_ij = (rdot_i, rdot_j) */
44: PetscBLASInt m; /* matrix dimension */
45: PetscBLASInt n; /* matrix dimension */
46: PetscBLASInt nrhs; /* the number of right hand sides */
47: PetscBLASInt lda; /* the padded matrix dimension */
48: PetscBLASInt ldb; /* the padded vector dimension */
49: PetscReal *s; /* the singular values */
50: PetscReal rcond; /* the exit condition */
51: PetscBLASInt rank; /* the effective rank */
52: PetscScalar *work; /* the work vector */
53: PetscReal *rwork; /* the real work vector used for complex */
54: PetscBLASInt lwork; /* the size of the work vector */
55: PetscBLASInt info; /* the output condition */
57: PetscBool setup_called; /* indicates whether SNESSetUp_NGMRES() has been called */
58: } SNES_NGMRES;
60: #define H(i, j) ngmres->h[i * ngmres->msize + j]
61: #define Q(i, j) ngmres->q[i * ngmres->msize + j]
63: /* private functions that are shared components of the methods */
64: PETSC_INTERN PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES, PetscInt, PetscInt, Vec, PetscReal, Vec);
65: PETSC_INTERN PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES, PetscInt, PetscInt, Vec, Vec, PetscReal, Vec, Vec, Vec);
66: PETSC_INTERN PetscErrorCode SNESNGMRESNorms_Private(SNES, PetscInt, Vec, Vec, Vec, Vec, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
67: PETSC_INTERN PetscErrorCode SNESNGMRESSelect_Private(SNES, PetscInt, Vec, Vec, PetscReal, PetscReal, PetscReal, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
68: PETSC_INTERN PetscErrorCode SNESNGMRESSelectRestart_Private(SNES, PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscBool *);
70: PETSC_INTERN PetscErrorCode SNESDestroy_NGMRES(SNES);
71: PETSC_INTERN PetscErrorCode SNESReset_NGMRES(SNES);
72: PETSC_INTERN PetscErrorCode SNESSetUp_NGMRES(SNES);
73: PETSC_INTERN PetscErrorCode SNESView_NGMRES(SNES, PetscViewer);