xref: /petsc/src/snes/impls/ngmres/snesngmres.h (revision 009bbdc485cd9ad46be9940d3549e2dde9cdc322)
1 #if !defined(_SNESNGMRES_H)
2 #define _SNESNGMRES_H
3 
4 #include <petsc-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   PetscInt    restart_periodic;/* number of iterations to restart after */
14 
15   SNESNGMRESRestartType   restart_type;
16   SNESNGMRESSelectType    select_type;
17 
18   /* History and subspace data */
19   Vec          *Fdot;          /* residual history -- length msize */
20   Vec          *Xdot;          /* solution history -- length msize */
21   PetscReal    *fnorms;        /* the residual norm history  */
22   PetscReal    *xnorms;        /* the solution norm history */
23 
24   /* General minimization problem context */
25   PetscScalar  *h;             /* the constraint matrix */
26   PetscScalar  *beta;          /* rhs for the minimization problem */
27   PetscScalar  *xi;            /* the dot-product of the current and previous res. */
28 
29   /* Line searches */
30   SNESLineSearch   additive_linesearch; /* Line search for the additive variant */
31 
32   /* Selection constants */
33   PetscBool    candidate;       /* use candidate storage approach */
34   PetscBool    singlereduction;/* use a single reduction (with more local work) for tolerance selection */
35   PetscReal    gammaA;         /* Criterion A residual tolerance */
36   PetscReal    epsilonB;       /* Criterion B difference tolerance */
37   PetscReal    deltaB;         /* Criterion B residual tolerance */
38   PetscReal    gammaC;         /* Restart residual tolerance */
39 
40   PetscReal    andersonBeta;   /* Relaxation parameter for Anderson Mixing */
41 
42   /* Least squares minimization solve context */
43   PetscScalar  *q;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
44   PetscBLASInt m;              /* matrix dimension */
45   PetscBLASInt n;              /* matrix dimension */
46   PetscBLASInt nrhs;           /* the number of right hand sides */
47   PetscBLASInt lda;            /* the padded matrix dimension */
48   PetscBLASInt ldb;            /* the padded vector dimension */
49   PetscReal    *s;             /* the singular values */
50   PetscReal    rcond;          /* the exit condition */
51   PetscBLASInt rank;           /* the effective rank */
52   PetscScalar  *work;          /* the work vector */
53   PetscReal    *rwork;         /* the real work vector used for complex */
54   PetscBLASInt lwork;          /* the size of the work vector */
55   PetscBLASInt info;           /* the output condition */
56 
57   PetscBool    setup_called;    /* indicates whether SNESSetUp_NGMRES() has been called  */
58 } SNES_NGMRES;
59 
60 #define H(i,j)  ngmres->h[i*ngmres->msize + j]
61 #define Q(i,j)  ngmres->q[i*ngmres->msize + j]
62 
63 /* private functions that are shared components of the methods */
64 PETSC_EXTERN PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES,PetscInt,PetscInt,Vec,PetscReal,Vec);
65 PETSC_EXTERN PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES,PetscInt,Vec,Vec,PetscReal,Vec,Vec,Vec);
66 PETSC_EXTERN PetscErrorCode SNESNGMRESCalculateDifferences_Private(SNES,PetscInt,Vec,Vec,Vec,Vec,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*);
67 PETSC_EXTERN PetscErrorCode SNESNGMRESSelect_Private(SNES,PetscInt,Vec,Vec,PetscReal,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscReal,Vec,Vec,Vec,PetscReal*);
68 PETSC_EXTERN PetscErrorCode SNESNGMRESSelectRestart_Private(SNES,PetscReal,PetscReal,PetscReal,PetscReal,PetscBool*);
69 
70 #endif
71