xref: /petsc/src/snes/impls/ngmres/snesngmres.h (revision fbada4e6d62d82be69a0720ca41ff3cebed5b071)
1 #ifndef _SNESNGMRES_H
2 #define _SNESNGMRES_H
3 
4 #include <private/snesimpl.h>
5 
6 /*  Data structure for the Nonlinear GMRES method.  */
7 typedef struct {
8 
9   /* Solver parameters and counters */
10   PetscInt     msize;          /* maximum size of krylov space */
11   PetscInt     k_rmax;         /* maximum number of iterations before restart */
12   PetscViewer  monitor;        /* debugging output for NGMRES */
13   /* History and subspace data */
14   Vec          *Fdot;          /* residual history -- length msize */
15   Vec          *Xdot;          /* solution history -- length msize */
16   PetscReal    *fnorms;        /* the residual norm history  */
17 
18   /* General minimization problem context */
19   PetscScalar  *h;             /* the constraint matrix */
20   PetscScalar  *beta;          /* rhs for the minimization problem */
21   PetscScalar  *xi;            /* the dot-product of the current and previous res. */
22 
23   /* Selection constants */
24   PetscReal    gammaA;         /* Criterion A residual tolerance */
25   PetscReal    epsilonB;       /* Criterion B difference tolerance */
26   PetscReal    deltaB;         /* Criterion B residual tolerance */
27   PetscReal    gammaC;         /* Restart residual tolerance */
28 
29   /* LS Minimization solve context */
30   PetscScalar  *q;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
31   PetscBLASInt m;              /* matrix dimension */
32   PetscBLASInt n;              /* matrix dimension */
33   PetscBLASInt nrhs;           /* the number of right hand sides */
34   PetscBLASInt lda;            /* the padded matrix dimension */
35   PetscBLASInt ldb;            /* the padded vector dimension */
36   PetscReal    *s;             /* the singular values */
37   PetscReal    rcond;          /* the exit condition */
38   PetscBLASInt rank;           /* the effective rank */
39   PetscScalar  *work;          /* the work vector */
40   PetscReal    *rwork;         /* the real work vector used for complex */
41   PetscBLASInt lwork;          /* the size of the work vector */
42   PetscBLASInt info;           /* the output condition */
43 
44 } SNES_NGMRES;
45 
46 #define H(i,j)  ngmres->h[i*ngmres->msize + j]
47 #define Q(i,j)  ngmres->q[i*ngmres->msize + j]
48 
49 #endif
50