xref: /petsc/src/snes/impls/ngmres/snesngmres.h (revision b549a33df52e18f44cd122c2139f7c5486754f82)
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   /* Selection constants */
25   PetscReal    gammaA;         /* Criterion A residual tolerance */
26   PetscReal    epsilonB;       /* Criterion B difference tolerance */
27   PetscReal    deltaB;         /* Criterion B residual tolerance */
28   PetscReal    gammaC;         /* Restart residual tolerance */
29 
30   /* LS Minimization solve context */
31   PetscScalar  *q;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
32   PetscBLASInt m;              /* matrix dimension */
33   PetscBLASInt n;              /* matrix dimension */
34   PetscBLASInt nrhs;           /* the number of right hand sides */
35   PetscBLASInt lda;            /* the padded matrix dimension */
36   PetscBLASInt ldb;            /* the padded vector dimension */
37   PetscReal    *s;             /* the singular values */
38   PetscReal    rcond;          /* the exit condition */
39   PetscBLASInt rank;           /* the effective rank */
40   PetscScalar  *work;          /* the work vector */
41   PetscReal    *rwork;         /* the real work vector used for complex */
42   PetscBLASInt lwork;          /* the size of the work vector */
43   PetscBLASInt info;           /* the output condition */
44 
45 } SNES_NGMRES;
46 
47 #define H(i,j)  ngmres->h[i*ngmres->msize + j]
48 #define Q(i,j)  ngmres->q[i*ngmres->msize + j]
49 
50 #endif
51