xref: /petsc/src/tao/bound/impls/bncg/bncg.h (revision 9dd11ecf0918283bb567d8b33a92f53ac4ea7840)
1*a4963045SJacob Faibussowitsch #pragma once
23ea99036SJacob Faibussowitsch 
3ac9112b8SAlp Dener /*
4ac9112b8SAlp Dener     Context for bound-constrained nonlinear conjugate gradient method
5ac9112b8SAlp Dener  */
6ac9112b8SAlp Dener 
7ac9112b8SAlp Dener #include <petsc/private/taoimpl.h>
8ac9112b8SAlp Dener 
9ac9112b8SAlp Dener typedef struct {
1050b47da0SAdam Denchfield   Mat B;
11484c7b14SAdam Denchfield   Mat pc;
1250b47da0SAdam Denchfield   Vec G_old, X_old, W, work;
13c8bcdf1eSAdam Denchfield   Vec g_work, y_work, d_work;
14c8bcdf1eSAdam Denchfield   Vec sk, yk;
1589da521bSAlp Dener   Vec unprojected_gradient, unprojected_gradient_old;
16ac9112b8SAlp Dener   Vec inactive_grad, inactive_step;
1761be54a6SAlp Dener 
18484c7b14SAdam Denchfield   IS active_lower, active_upper, active_fixed, active_idx, inactive_idx, inactive_old, new_inactives;
19484c7b14SAdam Denchfield 
2050b47da0SAdam Denchfield   PetscReal alpha; /* convex ratio in the scalar scaling */
2150b47da0SAdam Denchfield   PetscReal hz_theta;
2250b47da0SAdam Denchfield   PetscReal xi;    /* Parameter for KD, DK, and HZ methods. */
2350b47da0SAdam Denchfield   PetscReal theta; /* The convex combination parameter in the SSML_Broyden method. */
2450b47da0SAdam Denchfield   PetscReal zeta;  /* Another parameter, exclusive to Kou-Dai, modifying the y_k scalar contribution */
2550b47da0SAdam Denchfield   PetscReal hz_eta, dk_eta;
2650b47da0SAdam Denchfield   PetscReal bfgs_scale, dfp_scale; /* Scaling of the bfgs/dfp tau parameter in the bfgs and broyden methods. Default 1. */
2750b47da0SAdam Denchfield   PetscReal tau_bfgs, tau_dfp;
2850b47da0SAdam Denchfield   PetscReal as_step, as_tol, yts, yty, sts;
29484c7b14SAdam Denchfield   PetscReal f, delta_min, delta_max;
3050b47da0SAdam Denchfield   PetscReal epsilon; /* Machine precision unless changed by user */
31c8bcdf1eSAdam Denchfield   PetscReal eps_23;  /*  Two-thirds power of machine precision */
3250b47da0SAdam Denchfield 
33d6e07cdcSHong Zhang   TaoBNCGType cg_type;         /*  Formula to use */
3450b47da0SAdam Denchfield   PetscInt    min_restart_num; /* Restarts every x*n iterations, where n is the dimension */
35484c7b14SAdam Denchfield   PetscInt    ls_fails, resets, descent_error, skipped_updates, pure_gd_steps;
3650b47da0SAdam Denchfield   PetscInt    iter_quad, min_quad; /* Dynamic restart variables in Dai-Kou, SIAM J. Optim. Vol 23, pp. 296-320, Algorithm 4.1 */
3750b47da0SAdam Denchfield   PetscInt    as_type;
3850b47da0SAdam Denchfield 
39c8bcdf1eSAdam Denchfield   PetscBool inv_sig;
40c8bcdf1eSAdam Denchfield   PetscReal tol_quad;        /* tolerance for Dai-Kou dynamic restart */
41c8bcdf1eSAdam Denchfield   PetscBool dynamic_restart; /* Keeps track of whether or not to do a dynamic (KD) restart */
42c8bcdf1eSAdam Denchfield   PetscBool spaced_restart;  /* If true, restarts the CG method every x*n iterations */
43c8bcdf1eSAdam Denchfield   PetscBool use_dynamic_restart;
44c8bcdf1eSAdam Denchfield   PetscBool neg_xi;
4550b47da0SAdam Denchfield   PetscBool unscaled_restart; /* Gradient descent restarts are done without rescaling*/
46c8bcdf1eSAdam Denchfield   PetscBool diag_scaling;
47484c7b14SAdam Denchfield   PetscBool no_scaling;
4850b47da0SAdam Denchfield 
49ac9112b8SAlp Dener } TAO_BNCG;
50ac9112b8SAlp Dener 
5161be54a6SAlp Dener PETSC_INTERN PetscErrorCode TaoBNCGEstimateActiveSet(Tao, PetscInt);
52a1318120SAlp Dener PETSC_INTERN PetscErrorCode TaoBNCGBoundStep(Tao, PetscInt, Vec);
53c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal);
54c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao, PetscReal);
558ca2df50S PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscBool);
56c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGComputeDiagScaling(Tao, PetscReal, PetscReal);
57c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGResetUpdate(Tao, PetscReal);
58c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGCheckDynamicRestart(Tao, PetscReal, PetscReal, PetscReal, PetscBool *, PetscReal);
59