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