xref: /petsc/src/snes/impls/ngmres/snesngmres.h (revision ae9be28942d2d3eca8476af0b23adcfdc0f3733e)
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 
14   /* History and subspace data */
15   Vec          *Fdot;          /* residual history -- length msize */
16   Vec          *Xdot;          /* solution history -- length msize */
17   PetscReal    *fnorms;        /* the residual norm history  */
18 
19   /* General minimization problem context */
20   PetscScalar  *h;             /* the constraint matrix */
21   PetscScalar  *beta;          /* rhs for the minimization problem */
22   PetscScalar  *xi;            /* the dot-product of the current and previous res. */
23 
24   /* Line searches */
25   SNESLineSearch   additive_linesearch; /* Line search for the additive variant */
26 
27   /* Selection constants */
28   PetscBool    anderson;       /* use anderson-mixing approach */
29   PetscBool    additive;       /* use additive variant instead of selection */
30   PetscReal    gammaA;         /* Criterion A residual tolerance */
31   PetscReal    epsilonB;       /* Criterion B difference tolerance */
32   PetscReal    deltaB;         /* Criterion B residual tolerance */
33   PetscReal    gammaC;         /* Restart residual tolerance */
34 
35   /* Least squares minimization solve context */
36   PetscScalar  *q;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
37   PetscBLASInt m;              /* matrix dimension */
38   PetscBLASInt n;              /* matrix dimension */
39   PetscBLASInt nrhs;           /* the number of right hand sides */
40   PetscBLASInt lda;            /* the padded matrix dimension */
41   PetscBLASInt ldb;            /* the padded vector dimension */
42   PetscReal    *s;             /* the singular values */
43   PetscReal    rcond;          /* the exit condition */
44   PetscBLASInt rank;           /* the effective rank */
45   PetscScalar  *work;          /* the work vector */
46   PetscReal    *rwork;         /* the real work vector used for complex */
47   PetscBLASInt lwork;          /* the size of the work vector */
48   PetscBLASInt info;           /* the output condition */
49 
50   PetscBool    setup_called;    /* indicates whether SNESSetUp_NGMRES() has been called  */
51 } SNES_NGMRES;
52 
53 #define H(i,j)  ngmres->h[i*ngmres->msize + j]
54 #define Q(i,j)  ngmres->q[i*ngmres->msize + j]
55 
56 #endif
57