xref: /petsc/src/tao/bound/impls/bncg/bncg.h (revision 50b47da0f1ee44d44f7d7183a5d0083a5c160eaf)
1 /*
2     Context for bound-constrained nonlinear conjugate gradient method
3  */
4 
5 
6 #ifndef __TAO_BNCG_H
7 #define __TAO_BNCG_H
8 
9 #include <petsc/private/taoimpl.h>
10 
11 typedef struct {
12   Mat B;
13   Vec G_old, X_old, W, work;
14   Vec g_work, y_work, d_work;
15   Vec sk, yk;
16   Vec steffnm1, steffn, steffnp1;
17   Vec steffva, steffvatmp;
18   Vec unprojected_gradient, unprojected_gradient_old;
19   Vec inactive_grad, inactive_step;
20 
21   PetscReal alpha;                   /* convex ratio in the scalar scaling */
22   PetscReal hz_theta;
23   PetscReal xi;                      /* Parameter for KD, DK, and HZ methods. */
24   PetscReal theta;                   /* The convex combination parameter in the SSML_Broyden method. */
25   PetscReal zeta;                    /* Another parameter, exclusive to Kou-Dai, modifying the y_k scalar contribution */
26   PetscReal hz_eta, dk_eta;
27   PetscReal bfgs_scale, dfp_scale;   /* Scaling of the bfgs/dfp tau parameter in the bfgs and broyden methods. Default 1. */
28   PetscReal tau_bfgs, tau_dfp;
29   PetscReal as_step, as_tol, yts, yty, sts;
30   PetscReal f;
31   PetscReal epsilon;                 /* Machine precision unless changed by user */
32   PetscReal eps_23;                  /*  Two-thirds power of machine precision */
33 
34   IS  active_lower, active_upper, active_fixed, active_idx, inactive_idx, inactive_old, new_inactives;
35 
36   PetscInt cg_type;                  /*  Formula to use */
37   PetscInt  min_restart_num;         /* Restarts every x*n iterations, where n is the dimension */
38   PetscInt ls_fails, resets, descent_error;
39   PetscInt iter_quad, min_quad;      /* Dynamic restart variables in Dai-Kou, SIAM J. Optim. Vol 23, pp. 296-320, Algorithm 4.1 */
40   PetscInt  as_type;
41 
42   PetscBool recycle;
43   PetscBool inv_sig;
44   PetscReal tol_quad;                /* tolerance for Dai-Kou dynamic restart */
45   PetscBool dynamic_restart;         /* Keeps track of whether or not to do a dynamic (KD) restart */
46   PetscBool use_steffenson;          /* In development - attempting to use vector-based steffenson acceleration to the fixed point */
47   PetscBool spaced_restart;          /* If true, restarts the CG method every x*n iterations */
48   PetscBool use_dynamic_restart;
49   PetscBool neg_xi;
50   PetscBool unscaled_restart;        /* Gradient descent restarts are done without rescaling*/
51   PetscBool diag_scaling;
52 
53 } TAO_BNCG;
54 
55 #endif /* ifndef __TAO_BNCG_H */
56 
57 PETSC_INTERN PetscErrorCode TaoBNCGEstimateActiveSet(Tao, PetscInt);
58 PETSC_INTERN PetscErrorCode TaoBNCGBoundStep(Tao, PetscInt, Vec);
59 PETSC_EXTERN PetscErrorCode TaoBNCGSetRecycleFlag(Tao, PetscBool);
60 PETSC_INTERN PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal, PetscReal, PetscReal, PetscReal*, PetscReal);
61 PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao, PetscReal);
62 PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao, PetscReal, PetscReal, PetscReal, PetscReal, PetscBool, PetscReal, PetscReal);
63 PETSC_INTERN PetscErrorCode TaoBNCGComputeDiagScaling(Tao, PetscReal, PetscReal);
64 PETSC_INTERN PetscErrorCode TaoBNCGResetUpdate(Tao, PetscReal);
65 PETSC_INTERN PetscErrorCode TaoBNCGCheckDynamicRestart(Tao, PetscReal, PetscReal, PetscReal, PetscBool*, PetscReal);
66 PETSC_INTERN PetscErrorCode TaoBNCGSteffensonAcceleration(Tao);
67