xref: /petsc/src/tao/bound/impls/bncg/bncg.h (revision c8bcdf1ea017b917d33a7731f9bdf7d9fca0b1ea)
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   Vec G_old, X_old, U, V, W, work;
13   Vec g_work, y_work, d_work;
14   Vec sk, yk;
15   Vec steffnm1, steffn, steffnp1;
16   Vec steffva, steffvatmp, steffm1va, steffm2va;
17   Vec unprojected_gradient, unprojected_gradient_old;
18   IS  active_lower, active_upper, active_fixed, active_idx, inactive_idx, inactive_old, new_inactives;
19   Vec inactive_grad, inactive_step;
20   Vec diag; /* might be unused. Might have to use if Alp's approach with inverses isn't right */
21   Vec invD;
22   Vec invDnew;
23   PetscInt  as_type;
24   PetscReal as_step, as_tol, yts, yty, sts;
25 
26   PetscReal f;
27   PetscReal rho, pow;
28   PetscReal eta;         /*  Restart tolerance */
29   PetscReal epsilon;     /*  Machine precision  */
30   PetscReal eps_23;      /*  Two-thirds power of machine precision */
31   PetscBool recycle;
32   PetscBool inv_sig;
33   PetscInt cg_type;           /*  Formula to use */
34   PetscReal theta;        /* The convex combination parameter in the SSML_Broyden method. */
35   PetscReal gamma;          /* Stabilization parameter for multiple BNCG methods */
36   PetscReal zeta;         /* Another parameter, exclusive to Kou-Dai, modifying the y_k scalar contribution */
37   PetscReal bfgs_scale;   /* Scaling of the bfgs tau parameter in the bfgs and broyden methods. Default 1. */
38   PetscReal tau_bfgs, tau_dfp;
39 
40   Vec bfgs_work;
41   PetscReal dfp_scale;    /* Scaling of the dfp tau parameter in the dfp and broyden methods. Default 1. */
42   Vec dfp_work;
43 
44   PetscInt ls_fails, resets, broken_ortho, descent_error;
45   PetscInt iter_quad, min_quad;   /* Dynamic restart variables in Dai-Kou, SIAM J. Optim. Vol 23, pp. 296-320, Algorithm 4.1 */
46   PetscReal tol_quad;           /* tolerance for Dai-Kou dynamic restart */
47   PetscBool dynamic_restart;    /* Keeps track of whether or not to do a dynamic (KD) restart */
48   PetscBool use_steffenson;     /* In development - attempting to use vector-based steffenson acceleration to the fixed point */
49   PetscBool spaced_restart;     /* If true, restarts the CG method every x*n iterations */
50   PetscBool use_dynamic_restart;
51   PetscInt  min_restart_num;    /* Restarts every x*n iterations, where n is the dimension */
52   PetscBool neg_xi;
53   PetscBool unscaled_restart;
54   PetscBool diag_scaling;
55   PetscReal alpha; /* convex ratio in the scaling */
56   PetscReal beta; /* Exponential factor in the diagonal scaling */
57   PetscReal hz_theta;
58   PetscReal xi; /* Parameter for KD, DK, and HZ methods. */
59 } TAO_BNCG;
60 
61 #endif /* ifndef __TAO_BNCG_H */
62 
63 PETSC_INTERN PetscErrorCode TaoBNCGEstimateActiveSet(Tao, PetscInt);
64 PETSC_INTERN PetscErrorCode TaoBNCGBoundStep(Tao, PetscInt, Vec);
65 PETSC_EXTERN PetscErrorCode TaoBNCGSetRecycleFlag(Tao, PetscBool);
66 PETSC_INTERN PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal, PetscReal, PetscReal, PetscReal*, PetscReal);
67 PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao, PetscReal);
68 PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao, PetscReal, PetscReal, PetscReal, PetscReal, PetscBool, PetscReal, PetscReal);
69 PETSC_INTERN PetscErrorCode TaoBNCGComputeDiagScaling(Tao, PetscReal, PetscReal);
70 PETSC_INTERN PetscErrorCode TaoBNCGResetUpdate(Tao, PetscReal);
71 PETSC_INTERN PetscErrorCode TaoBNCGCheckDynamicRestart(Tao, PetscReal, PetscReal, PetscReal, PetscBool*, PetscReal);
72 PETSC_INTERN PetscErrorCode TaoBNCGSteffensonAcceleration(Tao);
73