Actual source code: nlsimpl.h
petsc-3.13.6 2020-09-29
1: /*
2: Context for a Newton line search method (unconstrained minimization)
3: */
7: #include <petsc/private/taoimpl.h>
9: typedef struct {
10: Mat M;
11: PC bfgs_pre;
13: Vec D;
14: Vec W;
16: Vec Xold;
17: Vec Gold;
19: /* Parameters when updating the perturbation added to the Hessian matrix
20: according to the following scheme:
22: pert = sval;
24: do until convergence
25: shift Hessian by pert
26: solve Newton system
28: if (linear solver failed or did not compute a descent direction)
29: use steepest descent direction and increase perturbation
31: if (0 == pert)
32: initialize perturbation
33: pert = min(imax, max(imin, imfac * norm(G)))
34: else
35: increase perturbation
36: pert = min(pmax, max(pgfac * pert, pmgfac * norm(G)))
37: fi
38: else
39: use linear solver direction and decrease perturbation
41: pert = min(psfac * pert, pmsfac * norm(G))
42: if (pert < pmin)
43: pert = 0
44: fi
45: fi
47: perform line search
48: function and gradient evaluation
49: check convergence
50: od
51: */
52: PetscReal sval; /* Starting perturbation value, default zero */
54: PetscReal imin; /* Minimum perturbation added during initialization */
55: PetscReal imax; /* Maximum perturbation added during initialization */
56: PetscReal imfac; /* Merit function factor during initialization */
58: PetscReal pmin; /* Minimim perturbation value */
59: PetscReal pmax; /* Maximum perturbation value */
60: PetscReal pgfac; /* Perturbation growth factor */
61: PetscReal psfac; /* Perturbation shrink factor */
62: PetscReal pmgfac; /* Merit function growth factor */
63: PetscReal pmsfac; /* Merit function shrink factor */
65: /* Parameters when updating the trust-region radius based on steplength
66: if step < nu1 (very bad step)
67: radius = omega1 * min(norm(d), radius)
68: elif step < nu2 (bad step)
69: radius = omega2 * min(norm(d), radius)
70: elif step < nu3 (okay step)
71: radius = omega3 * radius;
72: elif step < nu4 (good step)
73: radius = max(omega4 * norm(d), radius)
74: else (very good step)
75: radius = max(omega5 * norm(d), radius)
76: fi
77: */
78: PetscReal nu1; /* used to compute trust-region radius */
79: PetscReal nu2; /* used to compute trust-region radius */
80: PetscReal nu3; /* used to compute trust-region radius */
81: PetscReal nu4; /* used to compute trust-region radius */
83: PetscReal omega1; /* factor used for trust-region update */
84: PetscReal omega2; /* factor used for trust-region update */
85: PetscReal omega3; /* factor used for trust-region update */
86: PetscReal omega4; /* factor used for trust-region update */
87: PetscReal omega5; /* factor used for trust-region update */
89: /* Parameters when updating the trust-region radius based on reduction
91: kappa = ared / pred
92: if kappa < eta1 (very bad step)
93: radius = alpha1 * min(norm(d), radius)
94: elif kappa < eta2 (bad step)
95: radius = alpha2 * min(norm(d), radius)
96: elif kappa < eta3 (okay step)
97: radius = alpha3 * radius;
98: elif kappa < eta4 (good step)
99: radius = max(alpha4 * norm(d), radius)
100: else (very good step)
101: radius = max(alpha5 * norm(d), radius)
102: fi
103: */
104: PetscReal eta1; /* used to compute trust-region radius */
105: PetscReal eta2; /* used to compute trust-region radius */
106: PetscReal eta3; /* used to compute trust-region radius */
107: PetscReal eta4; /* used to compute trust-region radius */
109: PetscReal alpha1; /* factor used for trust-region update */
110: PetscReal alpha2; /* factor used for trust-region update */
111: PetscReal alpha3; /* factor used for trust-region update */
112: PetscReal alpha4; /* factor used for trust-region update */
113: PetscReal alpha5; /* factor used for trust-region update */
115: /* Parameters when updating the trust-region radius based on interpolation
117: kappa = ared / pred
118: if kappa >= 1.0 - mu1 (very good step)
119: choose tau in [gamma3, gamma4]
120: radius = max(tau * norm(d), radius)
121: elif kappa >= 1.0 - mu2 (good step)
122: choose tau in [gamma2, gamma3]
123: if (tau >= 1.0)
124: radius = max(tau * norm(d), radius)
125: else
126: radius = tau * min(norm(d), radius)
127: fi
128: else (bad step)
129: choose tau in [gamma1, 1.0]
130: radius = tau * min(norm(d), radius)
131: fi
132: */
133: PetscReal mu1; /* used for model agreement in interpolation */
134: PetscReal mu2; /* used for model agreement in interpolation */
136: PetscReal gamma1; /* factor used for interpolation */
137: PetscReal gamma2; /* factor used for interpolation */
138: PetscReal gamma3; /* factor used for interpolation */
139: PetscReal gamma4; /* factor used for interpolation */
141: PetscReal theta; /* factor used for interpolation */
143: /* Parameters when initializing trust-region radius based on interpolation */
144: PetscReal mu1_i; /* used for model agreement in interpolation */
145: PetscReal mu2_i; /* used for model agreement in interpolation */
147: PetscReal gamma1_i; /* factor used for interpolation */
148: PetscReal gamma2_i; /* factor used for interpolation */
149: PetscReal gamma3_i; /* factor used for interpolation */
150: PetscReal gamma4_i; /* factor used for interpolation */
152: PetscReal theta_i; /* factor used for interpolation */
154: /* Other parameters */
155: PetscReal min_radius; /* lower bound on initial radius value */
156: PetscReal max_radius; /* upper bound on trust region radius */
157: PetscReal epsilon; /* tolerance used when computing ared/pred */
159: PetscInt newt; /* Newton directions attempted */
160: PetscInt bfgs; /* BFGS directions attempted */
161: PetscInt grad; /* Gradient directions attempted */
163: PetscInt init_type; /* Trust-region initialization method */
164: PetscInt update_type; /* Trust-region update method */
166: PetscInt ksp_atol;
167: PetscInt ksp_rtol;
168: PetscInt ksp_ctol;
169: PetscInt ksp_negc;
170: PetscInt ksp_dtol;
171: PetscInt ksp_iter;
172: PetscInt ksp_othr;
173: } TAO_NLS;
175: #endif /* if !defined(__TAO_NLS_H) */