1fb90e4d1STodd Munson #include <../src/tao/unconstrained/impls/ntr/ntrimpl.h>
2a7e14dcfSSatish Balay
3aaa7dc30SBarry Smith #include <petscksp.h>
4a7e14dcfSSatish Balay
5a7e14dcfSSatish Balay #define NTR_INIT_CONSTANT 0
6a7e14dcfSSatish Balay #define NTR_INIT_DIRECTION 1
7a7e14dcfSSatish Balay #define NTR_INIT_INTERPOLATION 2
8a7e14dcfSSatish Balay #define NTR_INIT_TYPES 3
9a7e14dcfSSatish Balay
10a7e14dcfSSatish Balay #define NTR_UPDATE_REDUCTION 0
11a7e14dcfSSatish Balay #define NTR_UPDATE_INTERPOLATION 1
12a7e14dcfSSatish Balay #define NTR_UPDATE_TYPES 2
13a7e14dcfSSatish Balay
1453506e15SBarry Smith static const char *NTR_INIT[64] = {"constant", "direction", "interpolation"};
15a7e14dcfSSatish Balay
1653506e15SBarry Smith static const char *NTR_UPDATE[64] = {"reduction", "interpolation"};
17a7e14dcfSSatish Balay
18a7e14dcfSSatish Balay /*
19a7e14dcfSSatish Balay TaoSolve_NTR - Implements Newton's Method with a trust region approach
20a7e14dcfSSatish Balay for solving unconstrained minimization problems.
21a7e14dcfSSatish Balay
22a7e14dcfSSatish Balay The basic algorithm is taken from MINPACK-2 (dstrn).
23a7e14dcfSSatish Balay
24a7e14dcfSSatish Balay TaoSolve_NTR computes a local minimizer of a twice differentiable function
25a7e14dcfSSatish Balay f by applying a trust region variant of Newton's method. At each stage
26a7e14dcfSSatish Balay of the algorithm, we use the prconditioned conjugate gradient method to
27a7e14dcfSSatish Balay determine an approximate minimizer of the quadratic equation
28a7e14dcfSSatish Balay
29a7e14dcfSSatish Balay q(s) = <s, Hs + g>
30a7e14dcfSSatish Balay
31a7e14dcfSSatish Balay subject to the trust region constraint
32a7e14dcfSSatish Balay
33a7e14dcfSSatish Balay || s ||_M <= radius,
34a7e14dcfSSatish Balay
35a7e14dcfSSatish Balay where radius is the trust region radius and M is a symmetric positive
36a7e14dcfSSatish Balay definite matrix (the preconditioner). Here g is the gradient and H
37a7e14dcfSSatish Balay is the Hessian matrix.
38a7e14dcfSSatish Balay
3905de396fSBarry Smith Note: TaoSolve_NTR MUST use the iterative solver KSPNASH, KSPSTCG,
4005de396fSBarry Smith or KSPGLTR. Thus, we set KSPNASH, KSPSTCG, or KSPGLTR in this
41a7e14dcfSSatish Balay routine regardless of what the user may have previously specified.
42a7e14dcfSSatish Balay */
TaoSolve_NTR(Tao tao)43d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_NTR(Tao tao)
44d71ae5a4SJacob Faibussowitsch {
45a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data;
46fb90e4d1STodd Munson KSPType ksp_type;
470ad3a497SAlp Dener PetscBool is_nash, is_stcg, is_gltr, is_bfgs, is_jacobi, is_symmetric, sym_set;
48a7e14dcfSSatish Balay KSPConvergedReason ksp_reason;
49fb90e4d1STodd Munson PC pc;
50a7e14dcfSSatish Balay PetscReal fmin, ftrial, prered, actred, kappa, sigma, beta;
51a7e14dcfSSatish Balay PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
52a7e14dcfSSatish Balay PetscReal f, gnorm;
53a7e14dcfSSatish Balay
54a7e14dcfSSatish Balay PetscReal norm_d;
55a7e14dcfSSatish Balay PetscInt needH;
56a7e14dcfSSatish Balay
57a7e14dcfSSatish Balay PetscInt i_max = 5;
58a7e14dcfSSatish Balay PetscInt j_max = 1;
59a7e14dcfSSatish Balay PetscInt i, j, N, n, its;
60a7e14dcfSSatish Balay
61a7e14dcfSSatish Balay PetscFunctionBegin;
6248a46eb9SPierre Jolivet if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n"));
63a7e14dcfSSatish Balay
649566063dSJacob Faibussowitsch PetscCall(KSPGetType(tao->ksp, &ksp_type));
659566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ksp_type, KSPNASH, &is_nash));
669566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ksp_type, KSPSTCG, &is_stcg));
679566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ksp_type, KSPGLTR, &is_gltr));
683c859ba3SBarry Smith PetscCheck(is_nash || is_stcg || is_gltr, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "TAO_NTR requires nash, stcg, or gltr for the KSP");
69a7e14dcfSSatish Balay
70fb90e4d1STodd Munson /* Initialize the radius and modify if it is too large or small */
71fb90e4d1STodd Munson tao->trust = tao->trust0;
72a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius);
73a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius);
74a7e14dcfSSatish Balay
750c51296cSAlp Dener /* Allocate the vectors needed for the BFGS approximation */
769566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc));
779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
790c51296cSAlp Dener if (is_bfgs) {
800c51296cSAlp Dener tr->bfgs_pre = pc;
819566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(tr->bfgs_pre, &tr->M));
829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n));
839566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N));
849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(tr->M, n, n, N, N));
859566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(tr->M, tao->solution, tao->gradient));
869566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(tr->M, &sym_set, &is_symmetric));
873c859ba3SBarry Smith PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
881baa6e33SBarry Smith } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
89a7e14dcfSSatish Balay
90a7e14dcfSSatish Balay /* Check convergence criteria */
919566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
929566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
93*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
94a7e14dcfSSatish Balay needH = 1;
95a7e14dcfSSatish Balay
963ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
979566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
989566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0));
99dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1003ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
101a7e14dcfSSatish Balay
102a7e14dcfSSatish Balay /* Initialize trust-region radius */
103a7e14dcfSSatish Balay switch (tr->init_type) {
104a7e14dcfSSatish Balay case NTR_INIT_CONSTANT:
105a7e14dcfSSatish Balay /* Use the initial radius specified */
106a7e14dcfSSatish Balay break;
107a7e14dcfSSatish Balay
108a7e14dcfSSatish Balay case NTR_INIT_INTERPOLATION:
109a7e14dcfSSatish Balay /* Use the initial radius specified */
110a7e14dcfSSatish Balay max_radius = 0.0;
111a7e14dcfSSatish Balay
112a7e14dcfSSatish Balay for (j = 0; j < j_max; ++j) {
113a7e14dcfSSatish Balay fmin = f;
114a7e14dcfSSatish Balay sigma = 0.0;
115a7e14dcfSSatish Balay
116a7e14dcfSSatish Balay if (needH) {
1179566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
118a7e14dcfSSatish Balay needH = 0;
119a7e14dcfSSatish Balay }
120a7e14dcfSSatish Balay
121a7e14dcfSSatish Balay for (i = 0; i < i_max; ++i) {
1229566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tr->W));
1239566063dSJacob Faibussowitsch PetscCall(VecAXPY(tr->W, -tao->trust / gnorm, tao->gradient));
1249566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tr->W, &ftrial));
125a7e14dcfSSatish Balay
126a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) {
127a7e14dcfSSatish Balay tau = tr->gamma1_i;
1289371c9d4SSatish Balay } else {
129a7e14dcfSSatish Balay if (ftrial < fmin) {
130a7e14dcfSSatish Balay fmin = ftrial;
131a7e14dcfSSatish Balay sigma = -tao->trust / gnorm;
132a7e14dcfSSatish Balay }
133a7e14dcfSSatish Balay
1349566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->gradient, tao->stepdirection));
1359566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &prered));
136a7e14dcfSSatish Balay
137a7e14dcfSSatish Balay prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
138a7e14dcfSSatish Balay actred = f - ftrial;
1399371c9d4SSatish Balay if ((PetscAbsScalar(actred) <= tr->epsilon) && (PetscAbsScalar(prered) <= tr->epsilon)) {
140a7e14dcfSSatish Balay kappa = 1.0;
1419371c9d4SSatish Balay } else {
142a7e14dcfSSatish Balay kappa = actred / prered;
143a7e14dcfSSatish Balay }
144a7e14dcfSSatish Balay
145a7e14dcfSSatish Balay tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
146a7e14dcfSSatish Balay tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
147a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2);
148a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2);
149a7e14dcfSSatish Balay
15018cfbf8eSSatish Balay if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tr->mu1_i) {
151a7e14dcfSSatish Balay /* Great agreement */
152a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust);
153a7e14dcfSSatish Balay
154a7e14dcfSSatish Balay if (tau_max < 1.0) {
155a7e14dcfSSatish Balay tau = tr->gamma3_i;
1569371c9d4SSatish Balay } else if (tau_max > tr->gamma4_i) {
157a7e14dcfSSatish Balay tau = tr->gamma4_i;
1589371c9d4SSatish Balay } else {
159a7e14dcfSSatish Balay tau = tau_max;
160a7e14dcfSSatish Balay }
1619371c9d4SSatish Balay } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tr->mu2_i) {
162a7e14dcfSSatish Balay /* Good agreement */
163a7e14dcfSSatish Balay max_radius = PetscMax(max_radius, tao->trust);
164a7e14dcfSSatish Balay
165a7e14dcfSSatish Balay if (tau_max < tr->gamma2_i) {
166a7e14dcfSSatish Balay tau = tr->gamma2_i;
1679371c9d4SSatish Balay } else if (tau_max > tr->gamma3_i) {
168a7e14dcfSSatish Balay tau = tr->gamma3_i;
1699371c9d4SSatish Balay } else {
170a7e14dcfSSatish Balay tau = tau_max;
171a7e14dcfSSatish Balay }
1729371c9d4SSatish Balay } else {
173a7e14dcfSSatish Balay /* Not good agreement */
174a7e14dcfSSatish Balay if (tau_min > 1.0) {
175a7e14dcfSSatish Balay tau = tr->gamma2_i;
1769371c9d4SSatish Balay } else if (tau_max < tr->gamma1_i) {
177a7e14dcfSSatish Balay tau = tr->gamma1_i;
1789371c9d4SSatish Balay } else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
179a7e14dcfSSatish Balay tau = tr->gamma1_i;
1809371c9d4SSatish Balay } else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
181a7e14dcfSSatish Balay tau = tau_1;
1829371c9d4SSatish Balay } else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
183a7e14dcfSSatish Balay tau = tau_2;
1849371c9d4SSatish Balay } else {
185a7e14dcfSSatish Balay tau = tau_max;
186a7e14dcfSSatish Balay }
187a7e14dcfSSatish Balay }
188a7e14dcfSSatish Balay }
189a7e14dcfSSatish Balay tao->trust = tau * tao->trust;
190a7e14dcfSSatish Balay }
191a7e14dcfSSatish Balay
192a7e14dcfSSatish Balay if (fmin < f) {
193a7e14dcfSSatish Balay f = fmin;
1949566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
1959566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
196a7e14dcfSSatish Balay
1979566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
198*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
199a7e14dcfSSatish Balay needH = 1;
200a7e14dcfSSatish Balay
2019566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
2029566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, 1.0));
203dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
2043ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
205a7e14dcfSSatish Balay }
206a7e14dcfSSatish Balay }
207a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, max_radius);
208a7e14dcfSSatish Balay
209a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */
210a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius);
211a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius);
212a7e14dcfSSatish Balay break;
213a7e14dcfSSatish Balay
214a7e14dcfSSatish Balay default:
215a7e14dcfSSatish Balay /* Norm of the first direction will initialize radius */
216a7e14dcfSSatish Balay tao->trust = 0.0;
217a7e14dcfSSatish Balay break;
218a7e14dcfSSatish Balay }
219a7e14dcfSSatish Balay
220a7e14dcfSSatish Balay /* Have not converged; continue with Newton method */
2213ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) {
222e1e80dc8SAlp Dener /* Call general purpose update function */
223270bebe6SStefano Zampini if (tao->ops->update) {
224270bebe6SStefano Zampini PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
225270bebe6SStefano Zampini PetscCall(TaoComputeObjective(tao, tao->solution, &f));
226270bebe6SStefano Zampini }
2278931d482SJason Sarich ++tao->niter;
228ae93cb3cSJason Sarich tao->ksp_its = 0;
229a7e14dcfSSatish Balay /* Compute the Hessian */
230a7e14dcfSSatish Balay if (needH) {
2319566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
232a7e14dcfSSatish Balay needH = 0;
233a7e14dcfSSatish Balay }
234a7e14dcfSSatish Balay
2350c51296cSAlp Dener if (tr->bfgs_pre) {
236a7e14dcfSSatish Balay /* Update the limited memory preconditioner */
2379566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(tr->M, tao->solution, tao->gradient));
238a7e14dcfSSatish Balay }
239a7e14dcfSSatish Balay
2403ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) {
2419566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre));
242a7e14dcfSSatish Balay
243a7e14dcfSSatish Balay /* Solve the trust region subproblem */
2449566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
2459566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tao->gradient, tao->stepdirection));
2469566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its));
247a7e14dcfSSatish Balay tao->ksp_its += its;
248ae93cb3cSJason Sarich tao->ksp_tot_its += its;
2499566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
250a7e14dcfSSatish Balay
251a7e14dcfSSatish Balay if (0.0 == tao->trust) {
252a7e14dcfSSatish Balay /* Radius was uninitialized; use the norm of the direction */
253a7e14dcfSSatish Balay if (norm_d > 0.0) {
254a7e14dcfSSatish Balay tao->trust = norm_d;
255a7e14dcfSSatish Balay
256a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */
257a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius);
258a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius);
2599371c9d4SSatish Balay } else {
260a7e14dcfSSatish Balay /* The direction was bad; set radius to default value and re-solve
261a7e14dcfSSatish Balay the trust-region subproblem to get a direction */
262a7e14dcfSSatish Balay tao->trust = tao->trust0;
263a7e14dcfSSatish Balay
264a7e14dcfSSatish Balay /* Modify the radius if it is too large or small */
265a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->min_radius);
266a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius);
267a7e14dcfSSatish Balay
2689566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
2699566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tao->gradient, tao->stepdirection));
2709566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its));
271a7e14dcfSSatish Balay tao->ksp_its += its;
2722d9aa51bSJason Sarich tao->ksp_tot_its += its;
2739566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
274a7e14dcfSSatish Balay
2753c859ba3SBarry Smith PetscCheck(norm_d != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
276a7e14dcfSSatish Balay }
277a7e14dcfSSatish Balay }
2789566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
2799566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(tao->ksp, &ksp_reason));
2800c51296cSAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (tr->bfgs_pre)) {
281a7e14dcfSSatish Balay /* Preconditioner is numerically indefinite; reset the
282a7e14dcfSSatish Balay approximate if using BFGS preconditioning. */
2839566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(tr->M, PETSC_FALSE));
2849566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(tr->M, tao->solution, tao->gradient));
285a7e14dcfSSatish Balay }
286a7e14dcfSSatish Balay
287a7e14dcfSSatish Balay if (NTR_UPDATE_REDUCTION == tr->update_type) {
288a7e14dcfSSatish Balay /* Get predicted reduction */
2899566063dSJacob Faibussowitsch PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
290a7e14dcfSSatish Balay if (prered >= 0.0) {
291a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot
292a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should
293a7e14dcfSSatish Balay be rejected! */
294a7e14dcfSSatish Balay tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
2959371c9d4SSatish Balay } else {
296a7e14dcfSSatish Balay /* Compute trial step and function value */
2979566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tr->W));
2989566063dSJacob Faibussowitsch PetscCall(VecAXPY(tr->W, 1.0, tao->stepdirection));
2999566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tr->W, &ftrial));
300a7e14dcfSSatish Balay
301a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) {
302a7e14dcfSSatish Balay tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
303a7e14dcfSSatish Balay } else {
304a7e14dcfSSatish Balay /* Compute and actual reduction */
305a7e14dcfSSatish Balay actred = f - ftrial;
306a7e14dcfSSatish Balay prered = -prered;
3079371c9d4SSatish Balay if ((PetscAbsScalar(actred) <= tr->epsilon) && (PetscAbsScalar(prered) <= tr->epsilon)) {
308a7e14dcfSSatish Balay kappa = 1.0;
3099371c9d4SSatish Balay } else {
310a7e14dcfSSatish Balay kappa = actred / prered;
311a7e14dcfSSatish Balay }
312a7e14dcfSSatish Balay
313a7e14dcfSSatish Balay /* Accept or reject the step and update radius */
314a7e14dcfSSatish Balay if (kappa < tr->eta1) {
315a7e14dcfSSatish Balay /* Reject the step */
316a7e14dcfSSatish Balay tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
3179371c9d4SSatish Balay } else {
318a7e14dcfSSatish Balay /* Accept the step */
319a7e14dcfSSatish Balay if (kappa < tr->eta2) {
320a7e14dcfSSatish Balay /* Marginal bad step */
321a7e14dcfSSatish Balay tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
3229371c9d4SSatish Balay } else if (kappa < tr->eta3) {
323a7e14dcfSSatish Balay /* Reasonable step */
324a7e14dcfSSatish Balay tao->trust = tr->alpha3 * tao->trust;
3259371c9d4SSatish Balay } else if (kappa < tr->eta4) {
326a7e14dcfSSatish Balay /* Good step */
327a7e14dcfSSatish Balay tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
3289371c9d4SSatish Balay } else {
329a7e14dcfSSatish Balay /* Very good step */
330a7e14dcfSSatish Balay tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
331a7e14dcfSSatish Balay }
332a7e14dcfSSatish Balay break;
333a7e14dcfSSatish Balay }
334a7e14dcfSSatish Balay }
335a7e14dcfSSatish Balay }
3369371c9d4SSatish Balay } else {
337a7e14dcfSSatish Balay /* Get predicted reduction */
3389566063dSJacob Faibussowitsch PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
339a7e14dcfSSatish Balay if (prered >= 0.0) {
340a7e14dcfSSatish Balay /* The predicted reduction has the wrong sign. This cannot
341a7e14dcfSSatish Balay happen in infinite precision arithmetic. Step should
342a7e14dcfSSatish Balay be rejected! */
343a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
3449371c9d4SSatish Balay } else {
3459566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tr->W));
3469566063dSJacob Faibussowitsch PetscCall(VecAXPY(tr->W, 1.0, tao->stepdirection));
3479566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tr->W, &ftrial));
348a7e14dcfSSatish Balay if (PetscIsInfOrNanReal(ftrial)) {
349a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
3509371c9d4SSatish Balay } else {
3519566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &beta));
352a7e14dcfSSatish Balay actred = f - ftrial;
353a7e14dcfSSatish Balay prered = -prered;
3549371c9d4SSatish Balay if ((PetscAbsScalar(actred) <= tr->epsilon) && (PetscAbsScalar(prered) <= tr->epsilon)) {
355a7e14dcfSSatish Balay kappa = 1.0;
3569371c9d4SSatish Balay } else {
357a7e14dcfSSatish Balay kappa = actred / prered;
358a7e14dcfSSatish Balay }
359a7e14dcfSSatish Balay
360a7e14dcfSSatish Balay tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
361a7e14dcfSSatish Balay tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
362a7e14dcfSSatish Balay tau_min = PetscMin(tau_1, tau_2);
363a7e14dcfSSatish Balay tau_max = PetscMax(tau_1, tau_2);
364a7e14dcfSSatish Balay
365a7e14dcfSSatish Balay if (kappa >= 1.0 - tr->mu1) {
366a7e14dcfSSatish Balay /* Great agreement; accept step and update radius */
367a7e14dcfSSatish Balay if (tau_max < 1.0) {
368a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
3699371c9d4SSatish Balay } else if (tau_max > tr->gamma4) {
370a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
3719371c9d4SSatish Balay } else {
372a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d);
373a7e14dcfSSatish Balay }
374a7e14dcfSSatish Balay break;
3759371c9d4SSatish Balay } else if (kappa >= 1.0 - tr->mu2) {
376a7e14dcfSSatish Balay /* Good agreement */
377a7e14dcfSSatish Balay
378a7e14dcfSSatish Balay if (tau_max < tr->gamma2) {
379a7e14dcfSSatish Balay tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
3809371c9d4SSatish Balay } else if (tau_max > tr->gamma3) {
381a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
3829371c9d4SSatish Balay } else if (tau_max < 1.0) {
383a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d);
3849371c9d4SSatish Balay } else {
385a7e14dcfSSatish Balay tao->trust = PetscMax(tao->trust, tau_max * norm_d);
386a7e14dcfSSatish Balay }
387a7e14dcfSSatish Balay break;
3889371c9d4SSatish Balay } else {
389a7e14dcfSSatish Balay /* Not good agreement */
390a7e14dcfSSatish Balay if (tau_min > 1.0) {
391a7e14dcfSSatish Balay tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
3929371c9d4SSatish Balay } else if (tau_max < tr->gamma1) {
393a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
3949371c9d4SSatish Balay } else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
395a7e14dcfSSatish Balay tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
3969371c9d4SSatish Balay } else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) && ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
397a7e14dcfSSatish Balay tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
3989371c9d4SSatish Balay } else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) && ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
399a7e14dcfSSatish Balay tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
4009371c9d4SSatish Balay } else {
401a7e14dcfSSatish Balay tao->trust = tau_max * PetscMin(tao->trust, norm_d);
402a7e14dcfSSatish Balay }
403a7e14dcfSSatish Balay }
404a7e14dcfSSatish Balay }
405a7e14dcfSSatish Balay }
406a7e14dcfSSatish Balay }
407a7e14dcfSSatish Balay
408a7e14dcfSSatish Balay /* The step computed was not good and the radius was decreased.
409a7e14dcfSSatish Balay Monitor the radius to terminate. */
4109566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
4119566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust));
412dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
413a7e14dcfSSatish Balay }
414a7e14dcfSSatish Balay
415a7e14dcfSSatish Balay /* The radius may have been increased; modify if it is too large */
416a7e14dcfSSatish Balay tao->trust = PetscMin(tao->trust, tr->max_radius);
417a7e14dcfSSatish Balay
4183ecd9318SAlp Dener if (tao->reason == TAO_CONTINUE_ITERATING) {
4199566063dSJacob Faibussowitsch PetscCall(VecCopy(tr->W, tao->solution));
420a7e14dcfSSatish Balay f = ftrial;
4219566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
4229566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
423*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
424a7e14dcfSSatish Balay needH = 1;
4259566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
4269566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, tao->trust));
427dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
428a7e14dcfSSatish Balay }
429a7e14dcfSSatish Balay }
4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
431a7e14dcfSSatish Balay }
432a7e14dcfSSatish Balay
TaoSetUp_NTR(Tao tao)433d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_NTR(Tao tao)
434d71ae5a4SJacob Faibussowitsch {
435a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data;
436a7e14dcfSSatish Balay
437a7e14dcfSSatish Balay PetscFunctionBegin;
4389566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
4399566063dSJacob Faibussowitsch if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
4409566063dSJacob Faibussowitsch if (!tr->W) PetscCall(VecDuplicate(tao->solution, &tr->W));
441a7e14dcfSSatish Balay
44283c8fe1dSLisandro Dalcin tr->bfgs_pre = NULL;
44383c8fe1dSLisandro Dalcin tr->M = NULL;
4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
445a7e14dcfSSatish Balay }
446a7e14dcfSSatish Balay
TaoDestroy_NTR(Tao tao)447d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_NTR(Tao tao)
448d71ae5a4SJacob Faibussowitsch {
449a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data;
450a7e14dcfSSatish Balay
451a7e14dcfSSatish Balay PetscFunctionBegin;
45248a46eb9SPierre Jolivet if (tao->setupcalled) PetscCall(VecDestroy(&tr->W));
453a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp));
4549566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
456a7e14dcfSSatish Balay }
457a7e14dcfSSatish Balay
TaoSetFromOptions_NTR(Tao tao,PetscOptionItems PetscOptionsObject)458ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_NTR(Tao tao, PetscOptionItems PetscOptionsObject)
459d71ae5a4SJacob Faibussowitsch {
460a7e14dcfSSatish Balay TAO_NTR *tr = (TAO_NTR *)tao->data;
461a7e14dcfSSatish Balay
462a7e14dcfSSatish Balay PetscFunctionBegin;
463d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Newton trust region method for unconstrained optimization");
4649566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_ntr_init_type", "tao->trust initialization type", "", NTR_INIT, NTR_INIT_TYPES, NTR_INIT[tr->init_type], &tr->init_type, NULL));
4659566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_ntr_update_type", "radius update type", "", NTR_UPDATE, NTR_UPDATE_TYPES, NTR_UPDATE[tr->update_type], &tr->update_type, NULL));
4669566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1, NULL));
4679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2, NULL));
4689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3, NULL));
4699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4, NULL));
4709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1, NULL));
4719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2, NULL));
4729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3, NULL));
4739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4, NULL));
4749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5, NULL));
4759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1, NULL));
4769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2, NULL));
4779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1, NULL));
4789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2, NULL));
4799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3, NULL));
4809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4, NULL));
4819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta, NULL));
4829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i, NULL));
4839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i, NULL));
4849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i, NULL));
4859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i, NULL));
4869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i, NULL));
4879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i, NULL));
4889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i, NULL));
4899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius, NULL));
4909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius, NULL));
4919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon, NULL));
492d0609cedSBarry Smith PetscOptionsHeadEnd();
4939566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp));
4943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
495a7e14dcfSSatish Balay }
496a7e14dcfSSatish Balay
4971522df2eSJason Sarich /*MC
4981522df2eSJason Sarich TAONTR - Newton's method with trust region for unconstrained minimization.
4991522df2eSJason Sarich At each iteration, the Newton trust region method solves the system.
5001522df2eSJason Sarich NTR expects a KSP solver with a trust region radius.
5011522df2eSJason Sarich min_d .5 dT Hk d + gkT d, s.t. ||d|| < Delta_k
5021522df2eSJason Sarich
5031522df2eSJason Sarich Options Database Keys:
5049d0a60b2SAlp Dener + -tao_ntr_init_type - "constant","direction","interpolation"
5051522df2eSJason Sarich . -tao_ntr_update_type - "reduction","interpolation"
5061522df2eSJason Sarich . -tao_ntr_min_radius - lower bound on trust region radius
5071522df2eSJason Sarich . -tao_ntr_max_radius - upper bound on trust region radius
5081522df2eSJason Sarich . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
5091522df2eSJason Sarich . -tao_ntr_mu1_i - mu1 interpolation init factor
5101522df2eSJason Sarich . -tao_ntr_mu2_i - mu2 interpolation init factor
5111522df2eSJason Sarich . -tao_ntr_gamma1_i - gamma1 interpolation init factor
5121522df2eSJason Sarich . -tao_ntr_gamma2_i - gamma2 interpolation init factor
5131522df2eSJason Sarich . -tao_ntr_gamma3_i - gamma3 interpolation init factor
5141522df2eSJason Sarich . -tao_ntr_gamma4_i - gamma4 interpolation init factor
5158966356dSPierre Jolivet . -tao_ntr_theta_i - theta1 interpolation init factor
5161522df2eSJason Sarich . -tao_ntr_eta1 - eta1 reduction update factor
5171522df2eSJason Sarich . -tao_ntr_eta2 - eta2 reduction update factor
5181522df2eSJason Sarich . -tao_ntr_eta3 - eta3 reduction update factor
5191522df2eSJason Sarich . -tao_ntr_eta4 - eta4 reduction update factor
5201522df2eSJason Sarich . -tao_ntr_alpha1 - alpha1 reduction update factor
5211522df2eSJason Sarich . -tao_ntr_alpha2 - alpha2 reduction update factor
5221522df2eSJason Sarich . -tao_ntr_alpha3 - alpha3 reduction update factor
5231522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
5241522df2eSJason Sarich . -tao_ntr_alpha4 - alpha4 reduction update factor
5251522df2eSJason Sarich . -tao_ntr_mu1 - mu1 interpolation update
5261522df2eSJason Sarich . -tao_ntr_mu2 - mu2 interpolation update
5271522df2eSJason Sarich . -tao_ntr_gamma1 - gamma1 interpolcation update
5281522df2eSJason Sarich . -tao_ntr_gamma2 - gamma2 interpolcation update
5291522df2eSJason Sarich . -tao_ntr_gamma3 - gamma3 interpolcation update
5301522df2eSJason Sarich . -tao_ntr_gamma4 - gamma4 interpolation update
5311522df2eSJason Sarich - -tao_ntr_theta - theta interpolation update
5321522df2eSJason Sarich
5331eb8069cSJason Sarich Level: beginner
5341522df2eSJason Sarich M*/
5351522df2eSJason Sarich
TaoCreate_NTR(Tao tao)536d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
537d71ae5a4SJacob Faibussowitsch {
538a7e14dcfSSatish Balay TAO_NTR *tr;
539a7e14dcfSSatish Balay
540a7e14dcfSSatish Balay PetscFunctionBegin;
5414dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&tr));
542a7e14dcfSSatish Balay
543a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_NTR;
544a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_NTR;
545a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_NTR;
546a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_NTR;
547a7e14dcfSSatish Balay
5486552cf8aSJason Sarich /* Override default settings (unless already changed) */
549606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
550606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 50);
551606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, trust0, 100.0);
552a7e14dcfSSatish Balay tao->data = (void *)tr;
553a7e14dcfSSatish Balay
554a7e14dcfSSatish Balay /* Standard trust region update parameters */
555a7e14dcfSSatish Balay tr->eta1 = 1.0e-4;
556a7e14dcfSSatish Balay tr->eta2 = 0.25;
557a7e14dcfSSatish Balay tr->eta3 = 0.50;
558a7e14dcfSSatish Balay tr->eta4 = 0.90;
559a7e14dcfSSatish Balay
560a7e14dcfSSatish Balay tr->alpha1 = 0.25;
561a7e14dcfSSatish Balay tr->alpha2 = 0.50;
562a7e14dcfSSatish Balay tr->alpha3 = 1.00;
563a7e14dcfSSatish Balay tr->alpha4 = 2.00;
564a7e14dcfSSatish Balay tr->alpha5 = 4.00;
565a7e14dcfSSatish Balay
566a7e14dcfSSatish Balay /* Interpolation trust region update parameters */
567a7e14dcfSSatish Balay tr->mu1 = 0.10;
568a7e14dcfSSatish Balay tr->mu2 = 0.50;
569a7e14dcfSSatish Balay
570a7e14dcfSSatish Balay tr->gamma1 = 0.25;
571a7e14dcfSSatish Balay tr->gamma2 = 0.50;
572a7e14dcfSSatish Balay tr->gamma3 = 2.00;
573a7e14dcfSSatish Balay tr->gamma4 = 4.00;
574a7e14dcfSSatish Balay
575a7e14dcfSSatish Balay tr->theta = 0.05;
576a7e14dcfSSatish Balay
577fb90e4d1STodd Munson /* Interpolation parameters for initialization */
578fb90e4d1STodd Munson tr->mu1_i = 0.35;
579fb90e4d1STodd Munson tr->mu2_i = 0.50;
580fb90e4d1STodd Munson
581fb90e4d1STodd Munson tr->gamma1_i = 0.0625;
582fb90e4d1STodd Munson tr->gamma2_i = 0.50;
583fb90e4d1STodd Munson tr->gamma3_i = 2.00;
584fb90e4d1STodd Munson tr->gamma4_i = 5.00;
585fb90e4d1STodd Munson
586fb90e4d1STodd Munson tr->theta_i = 0.25;
587fb90e4d1STodd Munson
588a7e14dcfSSatish Balay tr->min_radius = 1.0e-10;
589a7e14dcfSSatish Balay tr->max_radius = 1.0e10;
590a7e14dcfSSatish Balay tr->epsilon = 1.0e-6;
591a7e14dcfSSatish Balay
592a7e14dcfSSatish Balay tr->init_type = NTR_INIT_INTERPOLATION;
593a7e14dcfSSatish Balay tr->update_type = NTR_UPDATE_REDUCTION;
594a7e14dcfSSatish Balay
595a7e14dcfSSatish Balay /* Set linear solver to default for trust region */
5969566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
5979566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
5989566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
5999566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_ntr_"));
6009566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp, KSPSTCG));
6013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
602a7e14dcfSSatish Balay }
603