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