xref: /petsc/src/snes/impls/ngmres/snesngmres.h (revision 8cdbd7576f676fd2f0cf4ce91c9fbf7f9b0d4574)
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     restart_it;     /* number of iterations the restart conditions persist before restart */
12   PetscViewer  monitor;        /* debugging output for NGMRES */
13   PetscBool    useGS;          /* use nonlinear GS or not */
14 
15   /* History and subspace data */
16   Vec          *Fdot;          /* residual history -- length msize */
17   Vec          *Xdot;          /* solution history -- length msize */
18   PetscReal    *fnorms;        /* the residual norm history  */
19 
20   /* General minimization problem context */
21   PetscScalar  *h;             /* the constraint matrix */
22   PetscScalar  *beta;          /* rhs for the minimization problem */
23   PetscScalar  *xi;            /* the dot-product of the current and previous res. */
24 
25   /* Selection constants */
26   PetscReal    gammaA;         /* Criterion A residual tolerance */
27   PetscReal    epsilonB;       /* Criterion B difference tolerance */
28   PetscReal    deltaB;         /* Criterion B residual tolerance */
29   PetscReal    gammaC;         /* Restart residual tolerance */
30 
31   /* LS Minimization solve context */
32   PetscScalar  *q;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
33   PetscBLASInt m;              /* matrix dimension */
34   PetscBLASInt n;              /* matrix dimension */
35   PetscBLASInt nrhs;           /* the number of right hand sides */
36   PetscBLASInt lda;            /* the padded matrix dimension */
37   PetscBLASInt ldb;            /* the padded vector dimension */
38   PetscReal    *s;             /* the singular values */
39   PetscReal    rcond;          /* the exit condition */
40   PetscBLASInt rank;           /* the effective rank */
41   PetscScalar  *work;          /* the work vector */
42   PetscReal    *rwork;         /* the real work vector used for complex */
43   PetscBLASInt lwork;          /* the size of the work vector */
44   PetscBLASInt info;           /* the output condition */
45 
46 } SNES_NGMRES;
47 
48 #define H(i,j)  ngmres->h[i*ngmres->msize + j]
49 #define Q(i,j)  ngmres->q[i*ngmres->msize + j]
50 
51 #endif
52