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