xref: /petsc/src/tao/unconstrained/impls/bmrm/bmrm.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1*a4963045SJacob Faibussowitsch #pragma once
2a7e14dcfSSatish Balay 
3af0996ceSBarry Smith #include <petsc/private/taoimpl.h>
4aaa7dc30SBarry Smith #include <petscmath.h>
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay #define BMRM_INFTY 1e30 /* single precision: ~\pm 10^{38.53}; PetscReal precision: ~\pm 10^{308.25} */
7a7e14dcfSSatish Balay #define ALPHA_MIN  1e-10
8a7e14dcfSSatish Balay #define ALPHA_MAX  1e10
9a7e14dcfSSatish Balay #define EPS_SV     1e-15
10a7e14dcfSSatish Balay #define EPS        1e-20
11a7e14dcfSSatish Balay #define TOL_LAM    1e-15
12a7e14dcfSSatish Balay #define TOL_R      1e-10
13a7e14dcfSSatish Balay #define INCRE_DIM  1000
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay /* Context for BMRM solver */
16a7e14dcfSSatish Balay typedef struct {
17a7e14dcfSSatish Balay   VecScatter scatter; /* Scatter context  */
18a7e14dcfSSatish Balay   Vec        local_w;
19a7e14dcfSSatish Balay   PetscReal  lambda;
20a7e14dcfSSatish Balay } TAO_BMRM;
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay typedef struct Vec_Chain {
23a7e14dcfSSatish Balay   Vec               V;
24a7e14dcfSSatish Balay   struct Vec_Chain *next;
25a7e14dcfSSatish Balay } Vec_Chain;
26a7e14dcfSSatish Balay 
27a7e14dcfSSatish Balay /* Context for Dai-Fletcher solver */
28a7e14dcfSSatish Balay typedef struct {
29a7e14dcfSSatish Balay   PetscInt   maxProjIter;
30a7e14dcfSSatish Balay   PetscInt   maxPGMIter;
31a7e14dcfSSatish Balay   PetscInt  *ipt, *ipt2, *uv;
32a7e14dcfSSatish Balay   PetscReal *g, *y, *tempv, *d, *Qd, *t, *xplus, *tplus, *sk, *yk;
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   PetscInt dim;
35a7e14dcfSSatish Balay 
36a7e14dcfSSatish Balay   PetscInt cur_num_cp;
37a7e14dcfSSatish Balay 
38a7e14dcfSSatish Balay   /* Variables (i.e. Lagrangian multipliers) */
39a7e14dcfSSatish Balay   PetscReal *x;
40a7e14dcfSSatish Balay 
41a7e14dcfSSatish Balay   /* Linear part of the objective function  */
42a7e14dcfSSatish Balay   PetscReal *f;
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   /* Hessian of the QP */
45a7e14dcfSSatish Balay   PetscReal **Q;
46a7e14dcfSSatish Balay 
47a7e14dcfSSatish Balay   /* Constraint matrix  */
48a7e14dcfSSatish Balay   PetscReal *a;
49a7e14dcfSSatish Balay 
50a7e14dcfSSatish Balay   /* RHS of the equality constraint */
51a7e14dcfSSatish Balay   PetscReal b;
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay   /* Lower bound vector for the variables */
54a7e14dcfSSatish Balay   PetscReal *l;
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay   /* Upper bound vector for the variables */
57a7e14dcfSSatish Balay   PetscReal *u;
58a7e14dcfSSatish Balay 
59a7e14dcfSSatish Balay   /* Tolerance for optimization error */
60a7e14dcfSSatish Balay   PetscReal tol;
61a7e14dcfSSatish Balay } TAO_DF;
62