xref: /petsc/src/snes/impls/ngmres/snesngmres.h (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
119653cdaSPeter Brune #ifndef _SNESNGMRES_H
219653cdaSPeter Brune #define _SNESNGMRES_H
319653cdaSPeter Brune 
4*b45d2f2cSJed Brown #include <petsc-private/snesimpl.h>
519653cdaSPeter Brune 
619653cdaSPeter Brune /*  Data structure for the Nonlinear GMRES method.  */
719653cdaSPeter Brune typedef struct {
819653cdaSPeter Brune 
998b3e84cSPeter Brune   /* Solver parameters and counters */
1098b3e84cSPeter Brune   PetscInt     msize;          /* maximum size of krylov space */
1128ed4a04SPeter Brune   PetscInt     restart_it;     /* number of iterations the restart conditions persist before restart */
12dfbf837cSBarry Smith   PetscViewer  monitor;        /* debugging output for NGMRES */
1328ed4a04SPeter Brune 
1498b3e84cSPeter Brune   /* History and subspace data */
15f109b39eSPeter Brune   Vec          *Fdot;          /* residual history -- length msize */
16f109b39eSPeter Brune   Vec          *Xdot;          /* solution history -- length msize */
17f109b39eSPeter Brune   PetscReal    *fnorms;        /* the residual norm history  */
1819653cdaSPeter Brune 
1998b3e84cSPeter Brune   /* General minimization problem context */
2098b3e84cSPeter Brune   PetscScalar  *h;             /* the constraint matrix */
2198b3e84cSPeter Brune   PetscScalar  *beta;          /* rhs for the minimization problem */
2298b3e84cSPeter Brune   PetscScalar  *xi;            /* the dot-product of the current and previous res. */
2319653cdaSPeter Brune 
24e7058c64SPeter Brune   /* Line searches */
25f1c6b773SPeter Brune   SNESLineSearch   additive_linesearch; /* Line search for the additive variant */
26e7058c64SPeter Brune 
2798b3e84cSPeter Brune   /* Selection constants */
28d2e16ddcSPeter Brune   PetscBool    anderson;       /* use anderson-mixing approach */
299f425c49SPeter Brune   PetscBool    additive;       /* use additive variant instead of selection */
3098b3e84cSPeter Brune   PetscReal    gammaA;         /* Criterion A residual tolerance */
3198b3e84cSPeter Brune   PetscReal    epsilonB;       /* Criterion B difference tolerance */
3298b3e84cSPeter Brune   PetscReal    deltaB;         /* Criterion B residual tolerance */
3398b3e84cSPeter Brune   PetscReal    gammaC;         /* Restart residual tolerance */
3419653cdaSPeter Brune 
35e7058c64SPeter Brune   /* Least squares minimization solve context */
3698b3e84cSPeter Brune   PetscScalar  *q;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
3798b3e84cSPeter Brune   PetscBLASInt m;              /* matrix dimension */
3898b3e84cSPeter Brune   PetscBLASInt n;              /* matrix dimension */
3998b3e84cSPeter Brune   PetscBLASInt nrhs;           /* the number of right hand sides */
4098b3e84cSPeter Brune   PetscBLASInt lda;            /* the padded matrix dimension */
4198b3e84cSPeter Brune   PetscBLASInt ldb;            /* the padded vector dimension */
4298b3e84cSPeter Brune   PetscReal    *s;             /* the singular values */
4398b3e84cSPeter Brune   PetscReal    rcond;          /* the exit condition */
4498b3e84cSPeter Brune   PetscBLASInt rank;           /* the effective rank */
4598b3e84cSPeter Brune   PetscScalar  *work;          /* the work vector */
4698b3e84cSPeter Brune   PetscReal    *rwork;         /* the real work vector used for complex */
4798b3e84cSPeter Brune   PetscBLASInt lwork;          /* the size of the work vector */
4898b3e84cSPeter Brune   PetscBLASInt info;           /* the output condition */
4919653cdaSPeter Brune 
5078440776SJed Brown   PetscBool    setup_called;    /* indicates whether SNESSetUp_NGMRES() has been called  */
5119653cdaSPeter Brune } SNES_NGMRES;
5219653cdaSPeter Brune 
5319653cdaSPeter Brune #define H(i,j)  ngmres->h[i*ngmres->msize + j]
5419653cdaSPeter Brune #define Q(i,j)  ngmres->q[i*ngmres->msize + j]
5519653cdaSPeter Brune 
5619653cdaSPeter Brune #endif
57