xref: /petsc/src/tao/unconstrained/impls/ntl/ntl.c (revision 03047865b8d8757cf1cf9cda45785c1537b01dc1)
155119615STodd Munson #include <../src/tao/unconstrained/impls/ntl/ntlimpl.h>
2a7e14dcfSSatish Balay 
3aaa7dc30SBarry Smith #include <petscksp.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay #define NTL_INIT_CONSTANT      0
6a7e14dcfSSatish Balay #define NTL_INIT_DIRECTION     1
7a7e14dcfSSatish Balay #define NTL_INIT_INTERPOLATION 2
8a7e14dcfSSatish Balay #define NTL_INIT_TYPES         3
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay #define NTL_UPDATE_REDUCTION     0
11a7e14dcfSSatish Balay #define NTL_UPDATE_INTERPOLATION 1
12a7e14dcfSSatish Balay #define NTL_UPDATE_TYPES         2
13a7e14dcfSSatish Balay 
1487f595a5SBarry Smith static const char *NTL_INIT[64] = {"constant", "direction", "interpolation"};
15a7e14dcfSSatish Balay 
1687f595a5SBarry Smith static const char *NTL_UPDATE[64] = {"reduction", "interpolation"};
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay /* Implements Newton's Method with a trust-region, line-search approach for
19a7e14dcfSSatish Balay    solving unconstrained minimization problems.  A More'-Thuente line search
20a7e14dcfSSatish Balay    is used to guarantee that the bfgs preconditioner remains positive
21a7e14dcfSSatish Balay    definite. */
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay #define NTL_NEWTON          0
24a7e14dcfSSatish Balay #define NTL_BFGS            1
25a7e14dcfSSatish Balay #define NTL_SCALED_GRADIENT 2
26a7e14dcfSSatish Balay #define NTL_GRADIENT        3
27a7e14dcfSSatish Balay 
TaoSolve_NTL(Tao tao)28d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_NTL(Tao tao)
29d71ae5a4SJacob Faibussowitsch {
30a7e14dcfSSatish Balay   TAO_NTL                     *tl = (TAO_NTL *)tao->data;
3155119615STodd Munson   KSPType                      ksp_type;
320ad3a497SAlp Dener   PetscBool                    is_nash, is_stcg, is_gltr, is_bfgs, is_jacobi, is_symmetric, sym_set;
33a7e14dcfSSatish Balay   KSPConvergedReason           ksp_reason;
3455119615STodd Munson   PC                           pc;
35e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay   PetscReal fmin, ftrial, prered, actred, kappa, sigma;
38a7e14dcfSSatish Balay   PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
39a7e14dcfSSatish Balay   PetscReal f, fold, gdx, gnorm;
40a7e14dcfSSatish Balay   PetscReal step = 1.0;
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay   PetscReal norm_d = 0.0;
43a7e14dcfSSatish Balay   PetscInt  stepType;
448931d482SJason Sarich   PetscInt  its;
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   PetscInt bfgsUpdates = 0;
47a7e14dcfSSatish Balay   PetscInt needH;
48a7e14dcfSSatish Balay 
49a7e14dcfSSatish Balay   PetscInt i_max = 5;
50a7e14dcfSSatish Balay   PetscInt j_max = 1;
51a7e14dcfSSatish Balay   PetscInt i, j, n, N;
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay   PetscInt tr_reject;
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay   PetscFunctionBegin;
5648a46eb9SPierre Jolivet   if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n"));
57a7e14dcfSSatish Balay 
589566063dSJacob Faibussowitsch   PetscCall(KSPGetType(tao->ksp, &ksp_type));
599566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPNASH, &is_nash));
609566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPSTCG, &is_stcg));
619566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPGLTR, &is_gltr));
623c859ba3SBarry Smith   PetscCheck(is_nash || is_stcg || is_gltr, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "TAO_NTR requires nash, stcg, or gltr for the KSP");
63a7e14dcfSSatish Balay 
6455119615STodd Munson   /* Initialize the radius and modify if it is too large or small */
6555119615STodd Munson   tao->trust = tao->trust0;
66a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tl->min_radius);
67a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tl->max_radius);
68a7e14dcfSSatish Balay 
690c51296cSAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
709566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
730c51296cSAlp Dener   if (is_bfgs) {
740c51296cSAlp Dener     tl->bfgs_pre = pc;
759566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(tl->bfgs_pre, &tl->M));
769566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
779566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
789566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(tl->M, n, n, N, N));
799566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(tl->M, tao->solution, tao->gradient));
809566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tl->M, &sym_set, &is_symmetric));
813c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
821baa6e33SBarry Smith   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay   /* Check convergence criteria */
859566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
869566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
87*76c63389SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
88a7e14dcfSSatish Balay   needH = 1;
89a7e14dcfSSatish Balay 
903ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
919566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
929566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
93dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
943ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
95a7e14dcfSSatish Balay 
9655119615STodd Munson   /* Initialize trust-region radius */
97a7e14dcfSSatish Balay   switch (tl->init_type) {
98a7e14dcfSSatish Balay   case NTL_INIT_CONSTANT:
99a7e14dcfSSatish Balay     /* Use the initial radius specified */
100a7e14dcfSSatish Balay     break;
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay   case NTL_INIT_INTERPOLATION:
103a7e14dcfSSatish Balay     /* Use the initial radius specified */
104a7e14dcfSSatish Balay     max_radius = 0.0;
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
107a7e14dcfSSatish Balay       fmin  = f;
108a7e14dcfSSatish Balay       sigma = 0.0;
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay       if (needH) {
1119566063dSJacob Faibussowitsch         PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
112a7e14dcfSSatish Balay         needH = 0;
113a7e14dcfSSatish Balay       }
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
1169566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->solution, tl->W));
1179566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tl->W, -tao->trust / gnorm, tao->gradient));
118a7e14dcfSSatish Balay 
1199566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, tl->W, &ftrial));
120a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
121a7e14dcfSSatish Balay           tau = tl->gamma1_i;
12253506e15SBarry Smith         } else {
123a7e14dcfSSatish Balay           if (ftrial < fmin) {
124a7e14dcfSSatish Balay             fmin  = ftrial;
125a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
126a7e14dcfSSatish Balay           }
127a7e14dcfSSatish Balay 
1289566063dSJacob Faibussowitsch           PetscCall(MatMult(tao->hessian, tao->gradient, tao->stepdirection));
1299566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &prered));
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
132a7e14dcfSSatish Balay           actred = f - ftrial;
13353506e15SBarry Smith           if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) {
134a7e14dcfSSatish Balay             kappa = 1.0;
13553506e15SBarry Smith           } else {
136a7e14dcfSSatish Balay             kappa = actred / prered;
137a7e14dcfSSatish Balay           }
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay           tau_1   = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred);
140a7e14dcfSSatish Balay           tau_2   = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred);
141a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
142a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
143a7e14dcfSSatish Balay 
14418cfbf8eSSatish Balay           if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tl->mu1_i) {
145a7e14dcfSSatish Balay             /* Great agreement */
146a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay             if (tau_max < 1.0) {
149a7e14dcfSSatish Balay               tau = tl->gamma3_i;
15053506e15SBarry Smith             } else if (tau_max > tl->gamma4_i) {
151a7e14dcfSSatish Balay               tau = tl->gamma4_i;
15253506e15SBarry Smith             } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) {
153a7e14dcfSSatish Balay               tau = tau_1;
15453506e15SBarry Smith             } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) {
155a7e14dcfSSatish Balay               tau = tau_2;
15653506e15SBarry Smith             } else {
157a7e14dcfSSatish Balay               tau = tau_max;
158a7e14dcfSSatish Balay             }
15918cfbf8eSSatish Balay           } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tl->mu2_i) {
160a7e14dcfSSatish Balay             /* Good agreement */
161a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
162a7e14dcfSSatish Balay 
163a7e14dcfSSatish Balay             if (tau_max < tl->gamma2_i) {
164a7e14dcfSSatish Balay               tau = tl->gamma2_i;
16553506e15SBarry Smith             } else if (tau_max > tl->gamma3_i) {
166a7e14dcfSSatish Balay               tau = tl->gamma3_i;
16753506e15SBarry Smith             } else {
168a7e14dcfSSatish Balay               tau = tau_max;
169a7e14dcfSSatish Balay             }
17053506e15SBarry Smith           } else {
171a7e14dcfSSatish Balay             /* Not good agreement */
172a7e14dcfSSatish Balay             if (tau_min > 1.0) {
173a7e14dcfSSatish Balay               tau = tl->gamma2_i;
17453506e15SBarry Smith             } else if (tau_max < tl->gamma1_i) {
175a7e14dcfSSatish Balay               tau = tl->gamma1_i;
17653506e15SBarry Smith             } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) {
177a7e14dcfSSatish Balay               tau = tl->gamma1_i;
17853506e15SBarry Smith             } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) {
179a7e14dcfSSatish Balay               tau = tau_1;
18053506e15SBarry Smith             } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) {
181a7e14dcfSSatish Balay               tau = tau_2;
18253506e15SBarry Smith             } else {
183a7e14dcfSSatish Balay               tau = tau_max;
184a7e14dcfSSatish Balay             }
185a7e14dcfSSatish Balay           }
186a7e14dcfSSatish Balay         }
187a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
188a7e14dcfSSatish Balay       }
189a7e14dcfSSatish Balay 
190a7e14dcfSSatish Balay       if (fmin < f) {
191a7e14dcfSSatish Balay         f = fmin;
1929566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
1939566063dSJacob Faibussowitsch         PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
194a7e14dcfSSatish Balay 
1959566063dSJacob Faibussowitsch         PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
196*76c63389SBarry Smith         PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
197a7e14dcfSSatish Balay         needH = 1;
198a7e14dcfSSatish Balay 
1999566063dSJacob Faibussowitsch         PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
2009566063dSJacob Faibussowitsch         PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
201dbbe0bcdSBarry Smith         PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
2023ba16761SJacob Faibussowitsch         if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
203a7e14dcfSSatish Balay       }
204a7e14dcfSSatish Balay     }
205a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay     /* Modify the radius if it is too large or small */
208a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tl->min_radius);
209a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tl->max_radius);
210a7e14dcfSSatish Balay     break;
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay   default:
213a7e14dcfSSatish Balay     /* Norm of the first direction will initialize radius */
214a7e14dcfSSatish Balay     tao->trust = 0.0;
215a7e14dcfSSatish Balay     break;
216a7e14dcfSSatish Balay   }
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
219a7e14dcfSSatish Balay   tl->ntrust = 0;
220a7e14dcfSSatish Balay   tl->newt   = 0;
221a7e14dcfSSatish Balay   tl->bfgs   = 0;
222a7e14dcfSSatish Balay   tl->grad   = 0;
223a7e14dcfSSatish Balay 
224a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2253ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
226e1e80dc8SAlp Dener     /* Call general purpose update function */
227270bebe6SStefano Zampini     if (tao->ops->update) {
228270bebe6SStefano Zampini       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
229270bebe6SStefano Zampini       PetscCall(TaoComputeObjective(tao, tao->solution, &f));
230270bebe6SStefano Zampini     }
2318931d482SJason Sarich     ++tao->niter;
232ae93cb3cSJason Sarich     tao->ksp_its = 0;
233a7e14dcfSSatish Balay     /* Compute the Hessian */
2341baa6e33SBarry Smith     if (needH) PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
235a7e14dcfSSatish Balay 
2360c51296cSAlp Dener     if (tl->bfgs_pre) {
237a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
2389566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
239a7e14dcfSSatish Balay       ++bfgsUpdates;
240a7e14dcfSSatish Balay     }
2419566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre));
242a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
2439566063dSJacob Faibussowitsch     PetscCall(KSPCGSetRadius(tao->ksp, tl->max_radius));
2449566063dSJacob Faibussowitsch     PetscCall(KSPSolve(tao->ksp, tao->gradient, tao->stepdirection));
2459566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp, &its));
246a7e14dcfSSatish Balay     tao->ksp_its += its;
247ae93cb3cSJason Sarich     tao->ksp_tot_its += its;
2489566063dSJacob Faibussowitsch     PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
249a7e14dcfSSatish Balay 
250a7e14dcfSSatish Balay     if (0.0 == tao->trust) {
251a7e14dcfSSatish Balay       /* Radius was uninitialized; use the norm of the direction */
252a7e14dcfSSatish Balay       if (norm_d > 0.0) {
253a7e14dcfSSatish Balay         tao->trust = norm_d;
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay         /* Modify the radius if it is too large or small */
256a7e14dcfSSatish Balay         tao->trust = PetscMax(tao->trust, tl->min_radius);
257a7e14dcfSSatish Balay         tao->trust = PetscMin(tao->trust, tl->max_radius);
25853506e15SBarry Smith       } else {
259a7e14dcfSSatish Balay         /* The direction was bad; set radius to default value and re-solve
260a7e14dcfSSatish Balay            the trust-region subproblem to get a direction */
261a7e14dcfSSatish Balay         tao->trust = tao->trust0;
262a7e14dcfSSatish Balay 
263a7e14dcfSSatish Balay         /* Modify the radius if it is too large or small */
264a7e14dcfSSatish Balay         tao->trust = PetscMax(tao->trust, tl->min_radius);
265a7e14dcfSSatish Balay         tao->trust = PetscMin(tao->trust, tl->max_radius);
266a7e14dcfSSatish Balay 
2679566063dSJacob Faibussowitsch         PetscCall(KSPCGSetRadius(tao->ksp, tl->max_radius));
2689566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, tao->gradient, tao->stepdirection));
2699566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp, &its));
270a7e14dcfSSatish Balay         tao->ksp_its += its;
2712d9aa51bSJason Sarich         tao->ksp_tot_its += its;
2729566063dSJacob Faibussowitsch         PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
273a7e14dcfSSatish Balay 
2743c859ba3SBarry Smith         PetscCheck(norm_d != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
275a7e14dcfSSatish Balay       }
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) && (tl->bfgs_pre)) {
281a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
282a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
2839566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
2849566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
285a7e14dcfSSatish Balay       bfgsUpdates = 1;
286a7e14dcfSSatish Balay     }
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay     /* Check trust-region reduction conditions */
289a7e14dcfSSatish Balay     tr_reject = 0;
290a7e14dcfSSatish Balay     if (NTL_UPDATE_REDUCTION == tl->update_type) {
291a7e14dcfSSatish Balay       /* Get predicted reduction */
2929566063dSJacob Faibussowitsch       PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
293a7e14dcfSSatish Balay       if (prered >= 0.0) {
294a7e14dcfSSatish Balay         /* The predicted reduction has the wrong sign.  This cannot
295a7e14dcfSSatish Balay            happen in infinite precision arithmetic.  Step should
296a7e14dcfSSatish Balay            be rejected! */
297a7e14dcfSSatish Balay         tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d);
298a7e14dcfSSatish Balay         tr_reject  = 1;
29953506e15SBarry Smith       } else {
300a7e14dcfSSatish Balay         /* Compute trial step and function value */
3019566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->solution, tl->W));
3029566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tl->W, 1.0, tao->stepdirection));
3039566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, tl->W, &ftrial));
304a7e14dcfSSatish Balay 
305a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
306a7e14dcfSSatish Balay           tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d);
307a7e14dcfSSatish Balay           tr_reject  = 1;
30853506e15SBarry Smith         } else {
309a7e14dcfSSatish Balay           /* Compute and actual reduction */
310a7e14dcfSSatish Balay           actred = f - ftrial;
311a7e14dcfSSatish Balay           prered = -prered;
3129371c9d4SSatish Balay           if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) {
313a7e14dcfSSatish Balay             kappa = 1.0;
31453506e15SBarry Smith           } else {
315a7e14dcfSSatish Balay             kappa = actred / prered;
316a7e14dcfSSatish Balay           }
317a7e14dcfSSatish Balay 
318a7e14dcfSSatish Balay           /* Accept of reject the step and update radius */
319a7e14dcfSSatish Balay           if (kappa < tl->eta1) {
320a7e14dcfSSatish Balay             /* Reject the step */
321a7e14dcfSSatish Balay             tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d);
322a7e14dcfSSatish Balay             tr_reject  = 1;
32353506e15SBarry Smith           } else {
324a7e14dcfSSatish Balay             /* Accept the step */
325a7e14dcfSSatish Balay             if (kappa < tl->eta2) {
326a7e14dcfSSatish Balay               /* Marginal bad step */
327a7e14dcfSSatish Balay               tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d);
32853506e15SBarry Smith             } else if (kappa < tl->eta3) {
329a7e14dcfSSatish Balay               /* Reasonable step */
330a7e14dcfSSatish Balay               tao->trust = tl->alpha3 * tao->trust;
33153506e15SBarry Smith             } else if (kappa < tl->eta4) {
332a7e14dcfSSatish Balay               /* Good step */
333a7e14dcfSSatish Balay               tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust);
33453506e15SBarry Smith             } else {
335a7e14dcfSSatish Balay               /* Very good step */
336a7e14dcfSSatish Balay               tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust);
337a7e14dcfSSatish Balay             }
338a7e14dcfSSatish Balay           }
339a7e14dcfSSatish Balay         }
340a7e14dcfSSatish Balay       }
34153506e15SBarry Smith     } else {
342a7e14dcfSSatish Balay       /* Get predicted reduction */
3439566063dSJacob Faibussowitsch       PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
344a7e14dcfSSatish Balay       if (prered >= 0.0) {
345a7e14dcfSSatish Balay         /* The predicted reduction has the wrong sign.  This cannot
346a7e14dcfSSatish Balay            happen in infinite precision arithmetic.  Step should
347a7e14dcfSSatish Balay            be rejected! */
348a7e14dcfSSatish Balay         tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
349a7e14dcfSSatish Balay         tr_reject  = 1;
35053506e15SBarry Smith       } else {
3519566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->solution, tl->W));
3529566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tl->W, 1.0, tao->stepdirection));
3539566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, tl->W, &ftrial));
354a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
355a7e14dcfSSatish Balay           tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
356a7e14dcfSSatish Balay           tr_reject  = 1;
35753506e15SBarry Smith         } else {
3589566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
359a7e14dcfSSatish Balay 
360a7e14dcfSSatish Balay           actred = f - ftrial;
361a7e14dcfSSatish Balay           prered = -prered;
3629371c9d4SSatish Balay           if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) {
363a7e14dcfSSatish Balay             kappa = 1.0;
36453506e15SBarry Smith           } else {
365a7e14dcfSSatish Balay             kappa = actred / prered;
366a7e14dcfSSatish Balay           }
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay           tau_1   = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred);
369a7e14dcfSSatish Balay           tau_2   = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred);
370a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
371a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
372a7e14dcfSSatish Balay 
373a7e14dcfSSatish Balay           if (kappa >= 1.0 - tl->mu1) {
374a7e14dcfSSatish Balay             /* Great agreement; accept step and update radius */
375a7e14dcfSSatish Balay             if (tau_max < 1.0) {
376a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d);
37753506e15SBarry Smith             } else if (tau_max > tl->gamma4) {
378a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d);
37953506e15SBarry Smith             } else {
380a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tau_max * norm_d);
381a7e14dcfSSatish Balay             }
38253506e15SBarry Smith           } else if (kappa >= 1.0 - tl->mu2) {
383a7e14dcfSSatish Balay             /* Good agreement */
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay             if (tau_max < tl->gamma2) {
386a7e14dcfSSatish Balay               tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d);
38753506e15SBarry Smith             } else if (tau_max > tl->gamma3) {
388a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d);
389a7e14dcfSSatish Balay             } else if (tau_max < 1.0) {
390a7e14dcfSSatish Balay               tao->trust = tau_max * PetscMin(tao->trust, norm_d);
39153506e15SBarry Smith             } else {
392a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tau_max * norm_d);
393a7e14dcfSSatish Balay             }
39453506e15SBarry Smith           } else {
395a7e14dcfSSatish Balay             /* Not good agreement */
396a7e14dcfSSatish Balay             if (tau_min > 1.0) {
397a7e14dcfSSatish Balay               tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d);
39853506e15SBarry Smith             } else if (tau_max < tl->gamma1) {
399a7e14dcfSSatish Balay               tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
40053506e15SBarry Smith             } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) {
401a7e14dcfSSatish Balay               tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
40253506e15SBarry Smith             } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) {
403a7e14dcfSSatish Balay               tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
40453506e15SBarry Smith             } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) {
405a7e14dcfSSatish Balay               tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
40653506e15SBarry Smith             } else {
407a7e14dcfSSatish Balay               tao->trust = tau_max * PetscMin(tao->trust, norm_d);
408a7e14dcfSSatish Balay             }
409a7e14dcfSSatish Balay             tr_reject = 1;
410a7e14dcfSSatish Balay           }
411a7e14dcfSSatish Balay         }
412a7e14dcfSSatish Balay       }
413a7e14dcfSSatish Balay     }
414a7e14dcfSSatish Balay 
415a7e14dcfSSatish Balay     if (tr_reject) {
416a7e14dcfSSatish Balay       /* The trust-region constraints rejected the step.  Apply a linesearch.
417a7e14dcfSSatish Balay          Check for descent direction. */
4189566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
419a7e14dcfSSatish Balay       if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
420*76c63389SBarry Smith         /* Newton step is not descent or direction produced infinity or NaN */
421a7e14dcfSSatish Balay 
4220c51296cSAlp Dener         if (!tl->bfgs_pre) {
423a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and updated
424a7e14dcfSSatish Balay              Must use gradient direction in this case */
4259566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->gradient, tao->stepdirection));
4269566063dSJacob Faibussowitsch           PetscCall(VecScale(tao->stepdirection, -1.0));
427a7e14dcfSSatish Balay           ++tl->grad;
428a7e14dcfSSatish Balay           stepType = NTL_GRADIENT;
42953506e15SBarry Smith         } else {
430a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
4319566063dSJacob Faibussowitsch           PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
4329566063dSJacob Faibussowitsch           PetscCall(VecScale(tao->stepdirection, -1.0));
433a7e14dcfSSatish Balay 
434a7e14dcfSSatish Balay           /* Check for success (descent direction) */
4359566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
436a7e14dcfSSatish Balay           if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
437a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
438a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case because
439a7e14dcfSSatish Balay                the first solve produces the scaled gradient direction,
440a7e14dcfSSatish Balay                which is guaranteed to be descent */
441a7e14dcfSSatish Balay 
442a7e14dcfSSatish Balay             /* Use steepest descent direction (scaled) */
4439566063dSJacob Faibussowitsch             PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
4449566063dSJacob Faibussowitsch             PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
4459566063dSJacob Faibussowitsch             PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
4469566063dSJacob Faibussowitsch             PetscCall(VecScale(tao->stepdirection, -1.0));
447a7e14dcfSSatish Balay 
448a7e14dcfSSatish Balay             bfgsUpdates = 1;
4490c51296cSAlp Dener             ++tl->grad;
4500c51296cSAlp Dener             stepType = NTL_GRADIENT;
45153506e15SBarry Smith           } else {
4529566063dSJacob Faibussowitsch             PetscCall(MatLMVMGetUpdateCount(tl->M, &bfgsUpdates));
453a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
454a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
4550c51296cSAlp Dener               ++tl->grad;
4560c51296cSAlp Dener               stepType = NTL_GRADIENT;
45753506e15SBarry Smith             } else {
458a7e14dcfSSatish Balay               ++tl->bfgs;
459a7e14dcfSSatish Balay               stepType = NTL_BFGS;
460a7e14dcfSSatish Balay             }
461a7e14dcfSSatish Balay           }
462a7e14dcfSSatish Balay         }
46353506e15SBarry Smith       } else {
464a7e14dcfSSatish Balay         /* Computed Newton step is descent */
465a7e14dcfSSatish Balay         ++tl->newt;
466a7e14dcfSSatish Balay         stepType = NTL_NEWTON;
467a7e14dcfSSatish Balay       }
468a7e14dcfSSatish Balay 
469a7e14dcfSSatish Balay       /* Perform the linesearch */
470a7e14dcfSSatish Balay       fold = f;
4719566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->solution, tl->Xold));
4729566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->gradient, tl->Gold));
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay       step = 1.0;
4759566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason));
4769566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
477a7e14dcfSSatish Balay 
47853506e15SBarry Smith       while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) { /* Linesearch failed */
479a7e14dcfSSatish Balay         /* Linesearch failed */
480a7e14dcfSSatish Balay         f = fold;
4819566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Xold, tao->solution));
4829566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Gold, tao->gradient));
483a7e14dcfSSatish Balay 
484a7e14dcfSSatish Balay         switch (stepType) {
485a7e14dcfSSatish Balay         case NTL_NEWTON:
486a7e14dcfSSatish Balay           /* Failed to obtain acceptable iterate with Newton step */
487a7e14dcfSSatish Balay 
4880c51296cSAlp Dener           if (tl->bfgs_pre) {
489a7e14dcfSSatish Balay             /* We don't have the bfgs matrix around and being updated
490a7e14dcfSSatish Balay                Must use gradient direction in this case */
4919566063dSJacob Faibussowitsch             PetscCall(VecCopy(tao->gradient, tao->stepdirection));
492a7e14dcfSSatish Balay             ++tl->grad;
493a7e14dcfSSatish Balay             stepType = NTL_GRADIENT;
49453506e15SBarry Smith           } else {
495a7e14dcfSSatish Balay             /* Attempt to use the BFGS direction */
4969566063dSJacob Faibussowitsch             PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
497a7e14dcfSSatish Balay 
498a7e14dcfSSatish Balay             /* Check for success (descent direction) */
4999566063dSJacob Faibussowitsch             PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
500a7e14dcfSSatish Balay             if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
501a7e14dcfSSatish Balay               /* BFGS direction is not descent or direction produced
502a7e14dcfSSatish Balay                  not a number.  We can assert bfgsUpdates > 1 in this case
503a7e14dcfSSatish Balay                  Use steepest descent direction (scaled) */
5049566063dSJacob Faibussowitsch               PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
5059566063dSJacob Faibussowitsch               PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
5069566063dSJacob Faibussowitsch               PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
507a7e14dcfSSatish Balay 
508a7e14dcfSSatish Balay               bfgsUpdates = 1;
5090c51296cSAlp Dener               ++tl->grad;
5100c51296cSAlp Dener               stepType = NTL_GRADIENT;
51153506e15SBarry Smith             } else {
5129566063dSJacob Faibussowitsch               PetscCall(MatLMVMGetUpdateCount(tl->M, &bfgsUpdates));
513a7e14dcfSSatish Balay               if (1 == bfgsUpdates) {
514a7e14dcfSSatish Balay                 /* The first BFGS direction is always the scaled gradient */
5150c51296cSAlp Dener                 ++tl->grad;
5160c51296cSAlp Dener                 stepType = NTL_GRADIENT;
51753506e15SBarry Smith               } else {
518a7e14dcfSSatish Balay                 ++tl->bfgs;
519a7e14dcfSSatish Balay                 stepType = NTL_BFGS;
520a7e14dcfSSatish Balay               }
521a7e14dcfSSatish Balay             }
522a7e14dcfSSatish Balay           }
523a7e14dcfSSatish Balay           break;
524a7e14dcfSSatish Balay 
525a7e14dcfSSatish Balay         case NTL_BFGS:
526a7e14dcfSSatish Balay           /* Can only enter if pc_type == NTL_PC_BFGS
527a7e14dcfSSatish Balay              Failed to obtain acceptable iterate with BFGS step
528a7e14dcfSSatish Balay              Attempt to use the scaled gradient direction */
5299566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
5309566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
5319566063dSJacob Faibussowitsch           PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
532a7e14dcfSSatish Balay 
533a7e14dcfSSatish Balay           bfgsUpdates = 1;
534a7e14dcfSSatish Balay           ++tl->grad;
535a7e14dcfSSatish Balay           stepType = NTL_GRADIENT;
536a7e14dcfSSatish Balay           break;
537a7e14dcfSSatish Balay         }
5389566063dSJacob Faibussowitsch         PetscCall(VecScale(tao->stepdirection, -1.0));
539a7e14dcfSSatish Balay 
540a7e14dcfSSatish Balay         /* This may be incorrect; linesearch has values for stepmax and stepmin
541a7e14dcfSSatish Balay            that should be reset. */
542a7e14dcfSSatish Balay         step = 1.0;
5439566063dSJacob Faibussowitsch         PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason));
5449566063dSJacob Faibussowitsch         PetscCall(TaoAddLineSearchCounts(tao));
545a7e14dcfSSatish Balay       }
546a7e14dcfSSatish Balay 
54753506e15SBarry Smith       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
548a7e14dcfSSatish Balay         /* Failed to find an improving point */
549a7e14dcfSSatish Balay         f = fold;
5509566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Xold, tao->solution));
5519566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Gold, tao->gradient));
552a7e14dcfSSatish Balay         tao->trust  = 0.0;
553a7e14dcfSSatish Balay         step        = 0.0;
554a7e14dcfSSatish Balay         tao->reason = TAO_DIVERGED_LS_FAILURE;
555a7e14dcfSSatish Balay         break;
55653506e15SBarry Smith       } else if (stepType == NTL_NEWTON) {
557a7e14dcfSSatish Balay         if (step < tl->nu1) {
558a7e14dcfSSatish Balay           /* Very bad step taken; reduce radius */
559a7e14dcfSSatish Balay           tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust);
56053506e15SBarry Smith         } else if (step < tl->nu2) {
561a7e14dcfSSatish Balay           /* Reasonably bad step taken; reduce radius */
562a7e14dcfSSatish Balay           tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust);
56353506e15SBarry Smith         } else if (step < tl->nu3) {
564a7e14dcfSSatish Balay           /* Reasonable step was taken; leave radius alone */
565a7e14dcfSSatish Balay           if (tl->omega3 < 1.0) {
566a7e14dcfSSatish Balay             tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust);
56753506e15SBarry Smith           } else if (tl->omega3 > 1.0) {
568a7e14dcfSSatish Balay             tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust);
569a7e14dcfSSatish Balay           }
57053506e15SBarry Smith         } else if (step < tl->nu4) {
571a7e14dcfSSatish Balay           /* Full step taken; increase the radius */
572a7e14dcfSSatish Balay           tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust);
57353506e15SBarry Smith         } else {
574a7e14dcfSSatish Balay           /* More than full step taken; increase the radius */
575a7e14dcfSSatish Balay           tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust);
576a7e14dcfSSatish Balay         }
57753506e15SBarry Smith       } else {
578a7e14dcfSSatish Balay         /* Newton step was not good; reduce the radius */
579a7e14dcfSSatish Balay         tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust);
580a7e14dcfSSatish Balay       }
58153506e15SBarry Smith     } else {
582a7e14dcfSSatish Balay       /* Trust-region step is accepted */
5839566063dSJacob Faibussowitsch       PetscCall(VecCopy(tl->W, tao->solution));
584a7e14dcfSSatish Balay       f = ftrial;
5859566063dSJacob Faibussowitsch       PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
586a7e14dcfSSatish Balay       ++tl->ntrust;
587a7e14dcfSSatish Balay     }
588a7e14dcfSSatish Balay 
589a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
590a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tl->max_radius);
591a7e14dcfSSatish Balay 
592e4cb33bbSBarry Smith     /* Check for converged */
5939566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
5943c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Not-a-Number");
595a7e14dcfSSatish Balay     needH = 1;
596a7e14dcfSSatish Balay 
5979566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
5989566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
599dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
600a7e14dcfSSatish Balay   }
6013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
602a7e14dcfSSatish Balay }
603a7e14dcfSSatish Balay 
TaoSetUp_NTL(Tao tao)604d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_NTL(Tao tao)
605d71ae5a4SJacob Faibussowitsch {
606a7e14dcfSSatish Balay   TAO_NTL *tl = (TAO_NTL *)tao->data;
607a7e14dcfSSatish Balay 
608a7e14dcfSSatish Balay   PetscFunctionBegin;
6099566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
6109566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
6119566063dSJacob Faibussowitsch   if (!tl->W) PetscCall(VecDuplicate(tao->solution, &tl->W));
6129566063dSJacob Faibussowitsch   if (!tl->Xold) PetscCall(VecDuplicate(tao->solution, &tl->Xold));
6139566063dSJacob Faibussowitsch   if (!tl->Gold) PetscCall(VecDuplicate(tao->solution, &tl->Gold));
61483c8fe1dSLisandro Dalcin   tl->bfgs_pre = NULL;
61583c8fe1dSLisandro Dalcin   tl->M        = NULL;
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617a7e14dcfSSatish Balay }
618a7e14dcfSSatish Balay 
TaoDestroy_NTL(Tao tao)619d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_NTL(Tao tao)
620d71ae5a4SJacob Faibussowitsch {
621a7e14dcfSSatish Balay   TAO_NTL *tl = (TAO_NTL *)tao->data;
622a7e14dcfSSatish Balay 
623a7e14dcfSSatish Balay   PetscFunctionBegin;
624a7e14dcfSSatish Balay   if (tao->setupcalled) {
6259566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tl->W));
6269566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tl->Xold));
6279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tl->Gold));
628a7e14dcfSSatish Balay   }
629a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
6309566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
6313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
632a7e14dcfSSatish Balay }
633a7e14dcfSSatish Balay 
TaoSetFromOptions_NTL(Tao tao,PetscOptionItems PetscOptionsObject)634ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_NTL(Tao tao, PetscOptionItems PetscOptionsObject)
635d71ae5a4SJacob Faibussowitsch {
636a7e14dcfSSatish Balay   TAO_NTL *tl = (TAO_NTL *)tao->data;
637a7e14dcfSSatish Balay 
638a7e14dcfSSatish Balay   PetscFunctionBegin;
639d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton trust region with line search method for unconstrained optimization");
6409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_ntl_init_type", "radius initialization type", "", NTL_INIT, NTL_INIT_TYPES, NTL_INIT[tl->init_type], &tl->init_type, NULL));
6419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_ntl_update_type", "radius update type", "", NTL_UPDATE, NTL_UPDATE_TYPES, NTL_UPDATE[tl->update_type], &tl->update_type, NULL));
6429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1, NULL));
6439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2, NULL));
6449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3, NULL));
6459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4, NULL));
6469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1, NULL));
6479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2, NULL));
6489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3, NULL));
6499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4, NULL));
6509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5, NULL));
6519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1, NULL));
6529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2, NULL));
6539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3, NULL));
6549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4, NULL));
6559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1, NULL));
6569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2, NULL));
6579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3, NULL));
6589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4, NULL));
6599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5, NULL));
6609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i, NULL));
6619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i, NULL));
6629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i, NULL));
6639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i, NULL));
6649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i, NULL));
6659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i, NULL));
6669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i, NULL));
6679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1, NULL));
6689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2, NULL));
6699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1, NULL));
6709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2, NULL));
6719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3, NULL));
6729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4, NULL));
6739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta, NULL));
6749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius, NULL));
6759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius, NULL));
6769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon, NULL));
677d0609cedSBarry Smith   PetscOptionsHeadEnd();
6789566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
6799566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
681a7e14dcfSSatish Balay }
682a7e14dcfSSatish Balay 
TaoView_NTL(Tao tao,PetscViewer viewer)683d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer)
684d71ae5a4SJacob Faibussowitsch {
685a7e14dcfSSatish Balay   TAO_NTL  *tl = (TAO_NTL *)tao->data;
686a7e14dcfSSatish Balay   PetscBool isascii;
687a7e14dcfSSatish Balay 
688a7e14dcfSSatish Balay   PetscFunctionBegin;
6899566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
690a7e14dcfSSatish Balay   if (isascii) {
6919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
69263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Trust-region steps: %" PetscInt_FMT "\n", tl->ntrust));
69363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton search steps: %" PetscInt_FMT "\n", tl->newt));
69463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS search steps: %" PetscInt_FMT "\n", tl->bfgs));
69563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient search steps: %" PetscInt_FMT "\n", tl->grad));
6969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
697a7e14dcfSSatish Balay   }
6983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
699a7e14dcfSSatish Balay }
700a7e14dcfSSatish Balay 
7011522df2eSJason Sarich /*MC
7023850be85SAlp Dener   TAONTL - Newton's method with trust region globalization and line search fallback.
7031522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system for d
7041522df2eSJason Sarich   and performs a line search in the d direction:
7051522df2eSJason Sarich 
7061522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
7071522df2eSJason Sarich 
7081522df2eSJason Sarich   Options Database Keys:
7099d0a60b2SAlp Dener + -tao_ntl_init_type - "constant","direction","interpolation"
7101522df2eSJason Sarich . -tao_ntl_update_type - "reduction","interpolation"
7111522df2eSJason Sarich . -tao_ntl_min_radius - lower bound on trust region radius
7121522df2eSJason Sarich . -tao_ntl_max_radius - upper bound on trust region radius
7131522df2eSJason Sarich . -tao_ntl_epsilon - tolerance for accepting actual / predicted reduction
7141522df2eSJason Sarich . -tao_ntl_mu1_i - mu1 interpolation init factor
7151522df2eSJason Sarich . -tao_ntl_mu2_i - mu2 interpolation init factor
7161522df2eSJason Sarich . -tao_ntl_gamma1_i - gamma1 interpolation init factor
7171522df2eSJason Sarich . -tao_ntl_gamma2_i - gamma2 interpolation init factor
7181522df2eSJason Sarich . -tao_ntl_gamma3_i - gamma3 interpolation init factor
7191522df2eSJason Sarich . -tao_ntl_gamma4_i - gamma4 interpolation init factor
7208966356dSPierre Jolivet . -tao_ntl_theta_i - theta1 interpolation init factor
7211522df2eSJason Sarich . -tao_ntl_eta1 - eta1 reduction update factor
7221522df2eSJason Sarich . -tao_ntl_eta2 - eta2 reduction update factor
7231522df2eSJason Sarich . -tao_ntl_eta3 - eta3 reduction update factor
7241522df2eSJason Sarich . -tao_ntl_eta4 - eta4 reduction update factor
7251522df2eSJason Sarich . -tao_ntl_alpha1 - alpha1 reduction update factor
7261522df2eSJason Sarich . -tao_ntl_alpha2 - alpha2 reduction update factor
7271522df2eSJason Sarich . -tao_ntl_alpha3 - alpha3 reduction update factor
7281522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor
7291522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor
7301522df2eSJason Sarich . -tao_ntl_mu1 - mu1 interpolation update
7311522df2eSJason Sarich . -tao_ntl_mu2 - mu2 interpolation update
7321522df2eSJason Sarich . -tao_ntl_gamma1 - gamma1 interpolcation update
7331522df2eSJason Sarich . -tao_ntl_gamma2 - gamma2 interpolcation update
7341522df2eSJason Sarich . -tao_ntl_gamma3 - gamma3 interpolcation update
7351522df2eSJason Sarich . -tao_ntl_gamma4 - gamma4 interpolation update
7361522df2eSJason Sarich - -tao_ntl_theta - theta1 interpolation update
7371522df2eSJason Sarich 
7381eb8069cSJason Sarich   Level: beginner
7391522df2eSJason Sarich M*/
TaoCreate_NTL(Tao tao)740d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_NTL(Tao tao)
741d71ae5a4SJacob Faibussowitsch {
742a7e14dcfSSatish Balay   TAO_NTL    *tl;
7438caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
74453506e15SBarry Smith 
745a7e14dcfSSatish Balay   PetscFunctionBegin;
7464dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&tl));
747a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_NTL;
748a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_NTL;
749a7e14dcfSSatish Balay   tao->ops->view           = TaoView_NTL;
750a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTL;
751a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_NTL;
752a7e14dcfSSatish Balay 
7536552cf8aSJason Sarich   /* Override default settings (unless already changed) */
754606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
755606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 50);
756606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, trust0, 100.0);
7576552cf8aSJason Sarich 
758a7e14dcfSSatish Balay   tao->data = (void *)tl;
759a7e14dcfSSatish Balay 
760a7e14dcfSSatish Balay   /* Default values for trust-region radius update based on steplength */
761a7e14dcfSSatish Balay   tl->nu1 = 0.25;
762a7e14dcfSSatish Balay   tl->nu2 = 0.50;
763a7e14dcfSSatish Balay   tl->nu3 = 1.00;
764a7e14dcfSSatish Balay   tl->nu4 = 1.25;
765a7e14dcfSSatish Balay 
766a7e14dcfSSatish Balay   tl->omega1 = 0.25;
767a7e14dcfSSatish Balay   tl->omega2 = 0.50;
768a7e14dcfSSatish Balay   tl->omega3 = 1.00;
769a7e14dcfSSatish Balay   tl->omega4 = 2.00;
770a7e14dcfSSatish Balay   tl->omega5 = 4.00;
771a7e14dcfSSatish Balay 
772a7e14dcfSSatish Balay   /* Default values for trust-region radius update based on reduction */
773a7e14dcfSSatish Balay   tl->eta1 = 1.0e-4;
774a7e14dcfSSatish Balay   tl->eta2 = 0.25;
775a7e14dcfSSatish Balay   tl->eta3 = 0.50;
776a7e14dcfSSatish Balay   tl->eta4 = 0.90;
777a7e14dcfSSatish Balay 
778a7e14dcfSSatish Balay   tl->alpha1 = 0.25;
779a7e14dcfSSatish Balay   tl->alpha2 = 0.50;
780a7e14dcfSSatish Balay   tl->alpha3 = 1.00;
781a7e14dcfSSatish Balay   tl->alpha4 = 2.00;
782a7e14dcfSSatish Balay   tl->alpha5 = 4.00;
783a7e14dcfSSatish Balay 
784a7e14dcfSSatish Balay   /* Default values for trust-region radius update based on interpolation */
785a7e14dcfSSatish Balay   tl->mu1 = 0.10;
786a7e14dcfSSatish Balay   tl->mu2 = 0.50;
787a7e14dcfSSatish Balay 
788a7e14dcfSSatish Balay   tl->gamma1 = 0.25;
789a7e14dcfSSatish Balay   tl->gamma2 = 0.50;
790a7e14dcfSSatish Balay   tl->gamma3 = 2.00;
791a7e14dcfSSatish Balay   tl->gamma4 = 4.00;
792a7e14dcfSSatish Balay 
793a7e14dcfSSatish Balay   tl->theta = 0.05;
794a7e14dcfSSatish Balay 
795a7e14dcfSSatish Balay   /* Default values for trust region initialization based on interpolation */
796a7e14dcfSSatish Balay   tl->mu1_i = 0.35;
797a7e14dcfSSatish Balay   tl->mu2_i = 0.50;
798a7e14dcfSSatish Balay 
799a7e14dcfSSatish Balay   tl->gamma1_i = 0.0625;
800a7e14dcfSSatish Balay   tl->gamma2_i = 0.5;
801a7e14dcfSSatish Balay   tl->gamma3_i = 2.0;
802a7e14dcfSSatish Balay   tl->gamma4_i = 5.0;
803a7e14dcfSSatish Balay 
804a7e14dcfSSatish Balay   tl->theta_i = 0.25;
805a7e14dcfSSatish Balay 
806a7e14dcfSSatish Balay   /* Remaining parameters */
807a7e14dcfSSatish Balay   tl->min_radius = 1.0e-10;
808a7e14dcfSSatish Balay   tl->max_radius = 1.0e10;
809a7e14dcfSSatish Balay   tl->epsilon    = 1.0e-6;
810a7e14dcfSSatish Balay 
811a7e14dcfSSatish Balay   tl->init_type   = NTL_INIT_INTERPOLATION;
812a7e14dcfSSatish Balay   tl->update_type = NTL_UPDATE_REDUCTION;
813a7e14dcfSSatish Balay 
8149566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
8159566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
8169566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
8179566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
8189566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
8199566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
8209566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
8219566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
8229566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_ntl_"));
8239566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
8243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
825a7e14dcfSSatish Balay }
826