Actual source code: ntrimpl.h
petsc-3.9.4 2018-09-11
1: /*
2: Context for a Newton trust region method (unconstrained minimization)
3: */
7: #include <petsc/private/taoimpl.h>
9: typedef struct {
10: Mat M;
12: Vec D;
13: Vec W;
15: Vec Diag;
16: PetscReal radius;
18: /* Parameters when updating the trust-region radius based on reduction
20: kappa = ared / pred
21: if kappa < eta1 (very bad step)
22: radius = alpha1 * min(norm(d), radius)
23: elif kappa < eta2 (bad step)
24: radius = alpha2 * min(norm(d), radius)
25: elif kappa < eta3 (okay step)
26: radius = alpha3 * radius;
27: elif kappa < eta4 (good step)
28: radius = max(alpha4 * norm(d), radius)
29: else (very good step)
30: radius = max(alpha5 * norm(d), radius)
31: fi
32: */
34: PetscReal eta1; /* used to compute trust-region radius */
35: PetscReal eta2; /* used to compute trust-region radius */
36: PetscReal eta3; /* used to compute trust-region radius */
37: PetscReal eta4; /* used to compute trust-region radius */
39: PetscReal alpha1; /* factor used for trust-region update */
40: PetscReal alpha2; /* factor used for trust-region update */
41: PetscReal alpha3; /* factor used for trust-region update */
42: PetscReal alpha4; /* factor used for trust-region update */
43: PetscReal alpha5; /* factor used for trust-region update */
45: /* Parameters when updating the trust-region radius based on interpolation
47: kappa = ared / pred
48: if kappa >= 1.0 - mu1 (very good step)
49: choose tau in [gamma3, gamma4]
50: radius = max(tau * norm(d), radius)
51: elif kappa >= 1.0 - mu2 (good step)
52: choose tau in [gamma2, gamma3]
53: if (tau >= 1.0)
54: radius = max(tau * norm(d), radius)
55: else
56: radius = tau * min(norm(d), radius)
57: fi
58: else (bad step)
59: choose tau in [gamma1, 1.0]
60: radius = tau * min(norm(d), radius)
61: fi
62: */
64: PetscReal mu1; /* used for model agreement in radius update */
65: PetscReal mu2; /* used for model agreement in radius update */
67: PetscReal gamma1; /* factor used for radius update */
68: PetscReal gamma2; /* factor used for radius update */
69: PetscReal gamma3; /* factor used for radius update */
70: PetscReal gamma4; /* factor used for radius update */
72: PetscReal theta; /* factor used for radius update */
74: /* Parameters when initializing trust-region radius based on interpolation */
76: PetscReal mu1_i; /* used for model agreement in interpolation */
77: PetscReal mu2_i; /* used for model agreement in interpolation */
79: PetscReal gamma1_i; /* factor used for interpolation */
80: PetscReal gamma2_i; /* factor used for interpolation */
81: PetscReal gamma3_i; /* factor used for interpolation */
82: PetscReal gamma4_i; /* factor used for interpolation */
84: PetscReal theta_i; /* factor used for interpolation */
86: PetscReal min_radius; /* lower bound on initial radius value */
87: PetscReal max_radius; /* upper bound on trust region radius */
88: PetscReal epsilon; /* tolerance used when computing actred/prered */
90: PetscInt pc_type; /* Preconditioner for the code */
91: PetscInt bfgs_scale_type;/* Scaling matrix for the bfgs preconditioner */
92: PetscInt init_type; /* Trust-region initialization method */
93: PetscInt update_type; /* Trust-region update method */
94: } TAO_NTR;
96: #endif /* if !defined(__TAO_NTR_H) */