Actual source code: nls.h
petsc-3.7.3 2016-08-01
1: /*
2: Context for a Newton line search method (unconstrained minimization)
3: */
5: #ifndef __TAO_NLS_H
7: #include <petsc/private/taoimpl.h>
9: typedef struct {
10: Mat M;
12: Vec D;
13: Vec W;
15: Vec Xold;
16: Vec Gold;
17: Vec Diag;
19: /* Parameters when updating the perturbation added to the Hessian matrix */
20: PetscReal sval; /* Starting perturbation value, default zero */
22: PetscReal imin; /* Minimum perturbation added during initialization */
23: PetscReal imax; /* Maximum perturbation added during initialization */
24: PetscReal imfac; /* Merit function factor during initialization */
26: PetscReal pmin; /* Minimim perturbation value */
27: PetscReal pmax; /* Maximum perturbation value */
28: PetscReal pgfac; /* Perturbation growth factor */
29: PetscReal psfac; /* Perturbation shrink factor */
30: PetscReal pmgfac; /* Merit function growth factor */
31: PetscReal pmsfac; /* Merit function shrink factor */
33: /* The perturbation to the Hessian matrix is initialized and updated
34: according to the following scheme:
36: pert = sval;
38: do until convergence
39: shift Hessian by pert
40: solve Newton system
42: if (linear solver failed or did not compute a descent direction)
43: use steepest descent direction and increase perturbation
45: if (0 == pert)
46: initialize perturbation
47: pert = min(imax, max(imin, imfac * norm(G)))
48: else
49: increase perturbation
50: pert = min(pmax, max(pgfac * pert, pmgfac * norm(G)))
51: fi
52: else
53: use linear solver direction and decrease perturbation
55: pert = min(psfac * pert, pmsfac * norm(G))
56: if (pert < pmin)
57: pert = 0
58: fi
59: fi
61: perform line search
62: function and gradient evaluation
63: check convergence
64: od
65: */
67: /* Parameters when updating the trust-region radius based on steplength */
68: PetscReal nu1; /* used to compute trust-region radius */
69: PetscReal nu2; /* used to compute trust-region radius */
70: PetscReal nu3; /* used to compute trust-region radius */
71: PetscReal nu4; /* used to compute trust-region radius */
73: PetscReal omega1; /* factor used for trust-region update */
74: PetscReal omega2; /* factor used for trust-region update */
75: PetscReal omega3; /* factor used for trust-region update */
76: PetscReal omega4; /* factor used for trust-region update */
77: PetscReal omega5; /* factor used for trust-region update */
79: /* if step < nu1 (very bad step)
80: radius = omega1 * min(norm(d), radius)
81: elif step < nu2 (bad step)
82: radius = omega2 * min(norm(d), radius)
83: elif step < nu3 (okay step)
84: radius = omega3 * radius;
85: elif step < nu4 (good step)
86: radius = max(omega4 * norm(d), radius)
87: else (very good step)
88: radius = max(omega5 * norm(d), radius)
89: fi
90: */
92: /* Parameters when updating the trust-region radius based on reduction */
93: PetscReal eta1; /* used to compute trust-region radius */
94: PetscReal eta2; /* used to compute trust-region radius */
95: PetscReal eta3; /* used to compute trust-region radius */
96: PetscReal eta4; /* used to compute trust-region radius */
98: PetscReal alpha1; /* factor used for trust-region update */
99: PetscReal alpha2; /* factor used for trust-region update */
100: PetscReal alpha3; /* factor used for trust-region update */
101: PetscReal alpha4; /* factor used for trust-region update */
102: PetscReal alpha5; /* factor used for trust-region update */
104: /* kappa = ared / pred
105: if kappa < eta1 (very bad step)
106: radius = alpha1 * min(norm(d), radius)
107: elif kappa < eta2 (bad step)
108: radius = alpha2 * min(norm(d), radius)
109: elif kappa < eta3 (okay step)
110: radius = alpha3 * radius;
111: elif kappa < eta4 (good step)
112: radius = max(alpha4 * norm(d), radius)
113: else (very good step)
114: radius = max(alpha5 * norm(d), radius)
115: fi
116: */
118: /* Parameters when updating the trust-region radius based on interpolation */
119: PetscReal mu1; /* used for model agreement in interpolation */
120: PetscReal mu2; /* used for model agreement in interpolation */
122: PetscReal gamma1; /* factor used for interpolation */
123: PetscReal gamma2; /* factor used for interpolation */
124: PetscReal gamma3; /* factor used for interpolation */
125: PetscReal gamma4; /* factor used for interpolation */
127: PetscReal theta; /* factor used for interpolation */
129: /* kappa = ared / pred
130: if kappa >= 1.0 - mu1 (very good step)
131: choose tau in [gamma3, gamma4]
132: radius = max(tau * norm(d), radius)
133: elif kappa >= 1.0 - mu2 (good step)
134: choose tau in [gamma2, gamma3]
135: if (tau >= 1.0)
136: radius = max(tau * norm(d), radius)
137: else
138: radius = tau * min(norm(d), radius)
139: fi
140: else (bad step)
141: choose tau in [gamma1, 1.0]
142: radius = tau * min(norm(d), radius)
143: fi
144: */
146: /* Parameters when initializing trust-region radius based on interpolation */
147: PetscReal mu1_i; /* used for model agreement in interpolation */
148: PetscReal mu2_i; /* used for model agreement in interpolation */
150: PetscReal gamma1_i; /* factor used for interpolation */
151: PetscReal gamma2_i; /* factor used for interpolation */
152: PetscReal gamma3_i; /* factor used for interpolation */
153: PetscReal gamma4_i; /* factor used for interpolation */
155: PetscReal theta_i; /* factor used for interpolation */
157: /* Other parameters */
158: PetscReal min_radius; /* lower bound on initial radius value */
159: PetscReal max_radius; /* upper bound on trust region radius */
160: PetscReal epsilon; /* tolerance used when computing ared/pred */
162: PetscInt newt; /* Newton directions attempted */
163: PetscInt bfgs; /* BFGS directions attempted */
164: PetscInt sgrad; /* Scaled gradient directions attempted */
165: PetscInt grad; /* Gradient directions attempted */
168: PetscInt ksp_type; /* KSP method for the code */
169: PetscInt pc_type; /* Preconditioner for the code */
170: PetscInt bfgs_scale_type; /* Scaling matrix to used for the bfgs preconditioner */
171: PetscInt init_type; /* Trust-region initialization method */
172: PetscInt update_type; /* Trust-region update method */
174: PetscInt ksp_atol;
175: PetscInt ksp_rtol;
176: PetscInt ksp_ctol;
177: PetscInt ksp_negc;
178: PetscInt ksp_dtol;
179: PetscInt ksp_iter;
180: PetscInt ksp_othr;
181: } TAO_NLS;
183: #endif /* ifndef __TAO_NLS_H */