xref: /petsc/src/tao/unconstrained/impls/ntl/ntlimpl.h (revision d2522c19e8fa9bca20aaca277941d9a63e71db6a)
1 /*
2   Context for a Newton trust-region, line-search method for unconstrained
3   minimization
4 */
5 
6 #if !defined(__TAO_NTL_H)
7 #define __TAO_NTL_H
8 #include <petsc/private/taoimpl.h>
9 
10 typedef struct {
11   Mat M;
12   PC  bfgs_pre;
13 
14   Vec W;
15   Vec Xold;
16   Vec Gold;
17 
18   /* Parameters when updating the trust-region radius based on steplength
19 
20      if   step < nu1            (very bad step)
21        radius = omega1 * min(norm(d), radius)
22      elif step < nu2            (bad step)
23        radius = omega2 * min(norm(d), radius)
24      elif step < nu3            (okay step)
25        radius = omega3 * radius;
26      elif step < nu4            (good step)
27        radius = max(omega4 * norm(d), radius)
28      else                       (very good step)
29        radius = max(omega5 * norm(d), radius)
30      fi
31   */
32 
33   PetscReal nu1; /* used to compute trust-region radius */
34   PetscReal nu2; /* used to compute trust-region radius */
35   PetscReal nu3; /* used to compute trust-region radius */
36   PetscReal nu4; /* used to compute trust-region radius */
37 
38   PetscReal omega1; /* factor used for trust-region update */
39   PetscReal omega2; /* factor used for trust-region update */
40   PetscReal omega3; /* factor used for trust-region update */
41   PetscReal omega4; /* factor used for trust-region update */
42   PetscReal omega5; /* factor used for trust-region update */
43 
44   /* Parameters when updating the trust-region radius based on reduction
45 
46      kappa = ared / pred
47      if   kappa < eta1          (very bad step)
48        radius = alpha1 * min(norm(d), radius)
49      elif kappa < eta2          (bad step)
50        radius = alpha2 * min(norm(d), radius)
51      elif kappa < eta3          (okay step)
52        radius = alpha3 * radius;
53      elif kappa < eta4          (good step)
54        radius = max(alpha4 * norm(d), radius)
55      else                       (very good step)
56        radius = max(alpha5 * norm(d), radius)
57      fi
58   */
59 
60   PetscReal eta1; /* used to compute trust-region radius */
61   PetscReal eta2; /* used to compute trust-region radius */
62   PetscReal eta3; /* used to compute trust-region radius */
63   PetscReal eta4; /* used to compute trust-region radius */
64 
65   PetscReal alpha1; /* factor used for trust-region update */
66   PetscReal alpha2; /* factor used for trust-region update */
67   PetscReal alpha3; /* factor used for trust-region update */
68   PetscReal alpha4; /* factor used for trust-region update */
69   PetscReal alpha5; /* factor used for trust-region update */
70 
71   /* Parameters when updating the trust-region radius based on interpolation
72      kappa = ared / pred
73      if   kappa >= 1.0 - mu1    (very good step)
74        choose tau in [gamma3, gamma4]
75        radius = max(tau * norm(d), radius)
76      elif kappa >= 1.0 - mu2    (good step)
77        choose tau in [gamma2, gamma3]
78        if (tau >= 1.0)
79          radius = max(tau * norm(d), radius)
80        else
81          radius = tau * min(norm(d), radius)
82        fi
83      else                       (bad step)
84        choose tau in [gamma1, 1.0]
85        radius = tau * min(norm(d), radius)
86      fi
87   */
88 
89   PetscReal mu1; /* used for model agreement in interpolation */
90   PetscReal mu2; /* used for model agreement in interpolation */
91 
92   PetscReal gamma1; /* factor used for interpolation */
93   PetscReal gamma2; /* factor used for interpolation */
94   PetscReal gamma3; /* factor used for interpolation */
95   PetscReal gamma4; /* factor used for interpolation */
96 
97   PetscReal theta; /* factor used for interpolation */
98 
99   /* Parameters when initializing trust-region radius based on interpolation */
100   PetscReal mu1_i; /* used for model agreement in interpolation */
101   PetscReal mu2_i; /* used for model agreement in interpolation */
102 
103   PetscReal gamma1_i; /* factor used for interpolation */
104   PetscReal gamma2_i; /* factor used for interpolation */
105   PetscReal gamma3_i; /* factor used for interpolation */
106   PetscReal gamma4_i; /* factor used for interpolation */
107 
108   PetscReal theta_i; /* factor used for interpolation */
109 
110   /* Other parameters */
111   PetscReal min_radius; /* lower bound on initial radius value */
112   PetscReal max_radius; /* upper bound on trust region radius */
113   PetscReal epsilon;    /* tolerance used when computing ared/pred */
114 
115   PetscInt ntrust; /* Trust-region steps accepted */
116   PetscInt newt;   /* Newton directions attempted */
117   PetscInt bfgs;   /* BFGS directions attempted */
118   PetscInt grad;   /* Gradient directions attempted */
119 
120   PetscInt init_type;   /* Trust-region initialization method */
121   PetscInt update_type; /* Trust-region update method */
122 } TAO_NTL;
123 
124 #endif /* if !defined(__TAO_NTL_H) */
125