xref: /petsc/src/tao/unconstrained/impls/ntl/ntl.c (revision a958fbfc1c07da5d8abfa22584ccb9c44e85e9ad)
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 
28441846f8SBarry Smith static PetscErrorCode TaoSolve_NTL(Tao tao)
29a7e14dcfSSatish Balay {
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;
56a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
579566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n"));
58a7e14dcfSSatish Balay   }
59a7e14dcfSSatish Balay 
609566063dSJacob Faibussowitsch   PetscCall(KSPGetType(tao->ksp,&ksp_type));
619566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type,KSPNASH,&is_nash));
629566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type,KSPSTCG,&is_stcg));
639566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type,KSPGLTR,&is_gltr));
643c859ba3SBarry Smith   PetscCheck(is_nash || is_stcg || is_gltr,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"TAO_NTR requires nash, stcg, or gltr for the KSP");
65a7e14dcfSSatish Balay 
6655119615STodd Munson   /* Initialize the radius and modify if it is too large or small */
6755119615STodd Munson   tao->trust = tao->trust0;
68a7e14dcfSSatish Balay   tao->trust = PetscMax(tao->trust, tl->min_radius);
69a7e14dcfSSatish Balay   tao->trust = PetscMin(tao->trust, tl->max_radius);
70a7e14dcfSSatish Balay 
710c51296cSAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
729566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
750c51296cSAlp Dener   if (is_bfgs) {
760c51296cSAlp Dener     tl->bfgs_pre = pc;
779566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(tl->bfgs_pre, &tl->M));
789566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
799566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
809566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(tl->M, n, n, N, N));
819566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(tl->M, tao->solution, tao->gradient));
829566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(tl->M, &sym_set, &is_symmetric));
833c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
840c51296cSAlp Dener   } else if (is_jacobi) {
859566063dSJacob Faibussowitsch     PetscCall(PCJacobiSetUseAbs(pc,PETSC_TRUE));
86a7e14dcfSSatish Balay   }
87a7e14dcfSSatish Balay 
88a7e14dcfSSatish Balay   /* Check convergence criteria */
899566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
909566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
913c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
92a7e14dcfSSatish Balay   needH = 1;
93a7e14dcfSSatish Balay 
943ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
959566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
969566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
979566063dSJacob Faibussowitsch   PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
983ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
99a7e14dcfSSatish Balay 
10055119615STodd Munson   /* Initialize trust-region radius */
101a7e14dcfSSatish Balay   switch(tl->init_type) {
102a7e14dcfSSatish Balay   case NTL_INIT_CONSTANT:
103a7e14dcfSSatish Balay     /* Use the initial radius specified */
104a7e14dcfSSatish Balay     break;
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay   case NTL_INIT_INTERPOLATION:
107a7e14dcfSSatish Balay     /* Use the initial radius specified */
108a7e14dcfSSatish Balay     max_radius = 0.0;
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay     for (j = 0; j < j_max; ++j) {
111a7e14dcfSSatish Balay       fmin = f;
112a7e14dcfSSatish Balay       sigma = 0.0;
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay       if (needH) {
1159566063dSJacob Faibussowitsch         PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre));
116a7e14dcfSSatish Balay         needH = 0;
117a7e14dcfSSatish Balay       }
118a7e14dcfSSatish Balay 
119a7e14dcfSSatish Balay       for (i = 0; i < i_max; ++i) {
1209566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->solution, tl->W));
1219566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tl->W, -tao->trust/gnorm, tao->gradient));
122a7e14dcfSSatish Balay 
1239566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, tl->W, &ftrial));
124a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
125a7e14dcfSSatish Balay           tau = tl->gamma1_i;
12653506e15SBarry Smith         } else {
127a7e14dcfSSatish Balay           if (ftrial < fmin) {
128a7e14dcfSSatish Balay             fmin = ftrial;
129a7e14dcfSSatish Balay             sigma = -tao->trust / gnorm;
130a7e14dcfSSatish Balay           }
131a7e14dcfSSatish Balay 
1329566063dSJacob Faibussowitsch           PetscCall(MatMult(tao->hessian, tao->gradient, tao->stepdirection));
1339566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &prered));
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
136a7e14dcfSSatish Balay           actred = f - ftrial;
13753506e15SBarry Smith           if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) {
138a7e14dcfSSatish Balay             kappa = 1.0;
13953506e15SBarry Smith           } else {
140a7e14dcfSSatish Balay             kappa = actred / prered;
141a7e14dcfSSatish Balay           }
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay           tau_1 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred);
144a7e14dcfSSatish Balay           tau_2 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred);
145a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
146a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
147a7e14dcfSSatish Balay 
14818cfbf8eSSatish Balay           if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tl->mu1_i) {
149a7e14dcfSSatish Balay             /* Great agreement */
150a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
151a7e14dcfSSatish Balay 
152a7e14dcfSSatish Balay             if (tau_max < 1.0) {
153a7e14dcfSSatish Balay               tau = tl->gamma3_i;
15453506e15SBarry Smith             } else if (tau_max > tl->gamma4_i) {
155a7e14dcfSSatish Balay               tau = tl->gamma4_i;
15653506e15SBarry Smith             } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) {
157a7e14dcfSSatish Balay               tau = tau_1;
15853506e15SBarry Smith             } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) {
159a7e14dcfSSatish Balay               tau = tau_2;
16053506e15SBarry Smith             } else {
161a7e14dcfSSatish Balay               tau = tau_max;
162a7e14dcfSSatish Balay             }
16318cfbf8eSSatish Balay           } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= tl->mu2_i) {
164a7e14dcfSSatish Balay             /* Good agreement */
165a7e14dcfSSatish Balay             max_radius = PetscMax(max_radius, tao->trust);
166a7e14dcfSSatish Balay 
167a7e14dcfSSatish Balay             if (tau_max < tl->gamma2_i) {
168a7e14dcfSSatish Balay               tau = tl->gamma2_i;
16953506e15SBarry Smith             } else if (tau_max > tl->gamma3_i) {
170a7e14dcfSSatish Balay               tau = tl->gamma3_i;
17153506e15SBarry Smith             } else {
172a7e14dcfSSatish Balay               tau = tau_max;
173a7e14dcfSSatish Balay             }
17453506e15SBarry Smith           } else {
175a7e14dcfSSatish Balay             /* Not good agreement */
176a7e14dcfSSatish Balay             if (tau_min > 1.0) {
177a7e14dcfSSatish Balay               tau = tl->gamma2_i;
17853506e15SBarry Smith             } else if (tau_max < tl->gamma1_i) {
179a7e14dcfSSatish Balay               tau = tl->gamma1_i;
18053506e15SBarry Smith             } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) {
181a7e14dcfSSatish Balay               tau = tl->gamma1_i;
18253506e15SBarry Smith             } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) &&  ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) {
183a7e14dcfSSatish Balay               tau = tau_1;
18453506e15SBarry Smith             } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) &&  ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) {
185a7e14dcfSSatish Balay               tau = tau_2;
18653506e15SBarry Smith             } else {
187a7e14dcfSSatish Balay               tau = tau_max;
188a7e14dcfSSatish Balay             }
189a7e14dcfSSatish Balay           }
190a7e14dcfSSatish Balay         }
191a7e14dcfSSatish Balay         tao->trust = tau * tao->trust;
192a7e14dcfSSatish Balay       }
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay       if (fmin < f) {
195a7e14dcfSSatish Balay         f = fmin;
1969566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
1979566063dSJacob Faibussowitsch         PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
198a7e14dcfSSatish Balay 
1999566063dSJacob Faibussowitsch         PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
2003c859ba3SBarry Smith         PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
201a7e14dcfSSatish Balay         needH = 1;
202a7e14dcfSSatish Balay 
2039566063dSJacob Faibussowitsch         PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
2049566063dSJacob Faibussowitsch         PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
2059566063dSJacob Faibussowitsch         PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
2063ecd9318SAlp Dener         if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
207a7e14dcfSSatish Balay       }
208a7e14dcfSSatish Balay     }
209a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, max_radius);
210a7e14dcfSSatish Balay 
211a7e14dcfSSatish Balay     /* Modify the radius if it is too large or small */
212a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, tl->min_radius);
213a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tl->max_radius);
214a7e14dcfSSatish Balay     break;
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay   default:
217a7e14dcfSSatish Balay     /* Norm of the first direction will initialize radius */
218a7e14dcfSSatish Balay     tao->trust = 0.0;
219a7e14dcfSSatish Balay     break;
220a7e14dcfSSatish Balay   }
221a7e14dcfSSatish Balay 
222a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
223a7e14dcfSSatish Balay   tl->ntrust = 0;
224a7e14dcfSSatish Balay   tl->newt = 0;
225a7e14dcfSSatish Balay   tl->bfgs = 0;
226a7e14dcfSSatish Balay   tl->grad = 0;
227a7e14dcfSSatish Balay 
228a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2293ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
230e1e80dc8SAlp Dener     /* Call general purpose update function */
231e1e80dc8SAlp Dener     if (tao->ops->update) {
2329566063dSJacob Faibussowitsch       PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update));
233e1e80dc8SAlp Dener     }
2348931d482SJason Sarich     ++tao->niter;
235ae93cb3cSJason Sarich     tao->ksp_its=0;
236a7e14dcfSSatish Balay     /* Compute the Hessian */
237a7e14dcfSSatish Balay     if (needH) {
2389566063dSJacob Faibussowitsch       PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre));
239a7e14dcfSSatish Balay     }
240a7e14dcfSSatish Balay 
2410c51296cSAlp Dener     if (tl->bfgs_pre) {
242a7e14dcfSSatish Balay       /* Update the limited memory preconditioner */
2439566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(tl->M,tao->solution, tao->gradient));
244a7e14dcfSSatish Balay       ++bfgsUpdates;
245a7e14dcfSSatish Balay     }
2469566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre));
247a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
2489566063dSJacob Faibussowitsch     PetscCall(KSPCGSetRadius(tao->ksp,tl->max_radius));
2499566063dSJacob Faibussowitsch     PetscCall(KSPSolve(tao->ksp, tao->gradient, tao->stepdirection));
2509566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(tao->ksp,&its));
251a7e14dcfSSatish Balay     tao->ksp_its+=its;
252ae93cb3cSJason Sarich     tao->ksp_tot_its+=its;
2539566063dSJacob Faibussowitsch     PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay     if (0.0 == tao->trust) {
256a7e14dcfSSatish Balay       /* Radius was uninitialized; use the norm of the direction */
257a7e14dcfSSatish Balay       if (norm_d > 0.0) {
258a7e14dcfSSatish Balay         tao->trust = norm_d;
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay         /* Modify the radius if it is too large or small */
261a7e14dcfSSatish Balay         tao->trust = PetscMax(tao->trust, tl->min_radius);
262a7e14dcfSSatish Balay         tao->trust = PetscMin(tao->trust, tl->max_radius);
26353506e15SBarry Smith       } else {
264a7e14dcfSSatish Balay         /* The direction was bad; set radius to default value and re-solve
265a7e14dcfSSatish Balay            the trust-region subproblem to get a direction */
266a7e14dcfSSatish Balay         tao->trust = tao->trust0;
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay         /* Modify the radius if it is too large or small */
269a7e14dcfSSatish Balay         tao->trust = PetscMax(tao->trust, tl->min_radius);
270a7e14dcfSSatish Balay         tao->trust = PetscMin(tao->trust, tl->max_radius);
271a7e14dcfSSatish Balay 
2729566063dSJacob Faibussowitsch         PetscCall(KSPCGSetRadius(tao->ksp,tl->max_radius));
2739566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, tao->gradient, tao->stepdirection));
2749566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp,&its));
275a7e14dcfSSatish Balay         tao->ksp_its+=its;
2762d9aa51bSJason Sarich         tao->ksp_tot_its+=its;
2779566063dSJacob Faibussowitsch         PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
278a7e14dcfSSatish Balay 
2793c859ba3SBarry Smith         PetscCheck(norm_d != 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero");
280a7e14dcfSSatish Balay       }
281a7e14dcfSSatish Balay     }
282a7e14dcfSSatish Balay 
2839566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
2849566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(tao->ksp, &ksp_reason));
2850c51296cSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (tl->bfgs_pre)) {
286a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
287a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
2889566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
2899566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
290a7e14dcfSSatish Balay       bfgsUpdates = 1;
291a7e14dcfSSatish Balay     }
292a7e14dcfSSatish Balay 
293a7e14dcfSSatish Balay     /* Check trust-region reduction conditions */
294a7e14dcfSSatish Balay     tr_reject = 0;
295a7e14dcfSSatish Balay     if (NTL_UPDATE_REDUCTION == tl->update_type) {
296a7e14dcfSSatish Balay       /* Get predicted reduction */
2979566063dSJacob Faibussowitsch       PetscCall(KSPCGGetObjFcn(tao->ksp,&prered));
298a7e14dcfSSatish Balay       if (prered >= 0.0) {
299a7e14dcfSSatish Balay         /* The predicted reduction has the wrong sign.  This cannot
300a7e14dcfSSatish Balay            happen in infinite precision arithmetic.  Step should
301a7e14dcfSSatish Balay            be rejected! */
302a7e14dcfSSatish Balay         tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d);
303a7e14dcfSSatish Balay         tr_reject = 1;
30453506e15SBarry Smith       } else {
305a7e14dcfSSatish Balay         /* Compute trial step and function value */
3069566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->solution, tl->W));
3079566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tl->W, 1.0, tao->stepdirection));
3089566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, tl->W, &ftrial));
309a7e14dcfSSatish Balay 
310a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
311a7e14dcfSSatish Balay           tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d);
312a7e14dcfSSatish Balay           tr_reject = 1;
31353506e15SBarry Smith         } else {
314a7e14dcfSSatish Balay           /* Compute and actual reduction */
315a7e14dcfSSatish Balay           actred = f - ftrial;
316a7e14dcfSSatish Balay           prered = -prered;
317a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tl->epsilon) &&
318a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tl->epsilon)) {
319a7e14dcfSSatish Balay             kappa = 1.0;
32053506e15SBarry Smith           } else {
321a7e14dcfSSatish Balay             kappa = actred / prered;
322a7e14dcfSSatish Balay           }
323a7e14dcfSSatish Balay 
324a7e14dcfSSatish Balay           /* Accept of reject the step and update radius */
325a7e14dcfSSatish Balay           if (kappa < tl->eta1) {
326a7e14dcfSSatish Balay             /* Reject the step */
327a7e14dcfSSatish Balay             tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d);
328a7e14dcfSSatish Balay             tr_reject = 1;
32953506e15SBarry Smith           } else {
330a7e14dcfSSatish Balay             /* Accept the step */
331a7e14dcfSSatish Balay             if (kappa < tl->eta2) {
332a7e14dcfSSatish Balay               /* Marginal bad step */
333a7e14dcfSSatish Balay               tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d);
33453506e15SBarry Smith             } else if (kappa < tl->eta3) {
335a7e14dcfSSatish Balay               /* Reasonable step */
336a7e14dcfSSatish Balay               tao->trust = tl->alpha3 * tao->trust;
33753506e15SBarry Smith             } else if (kappa < tl->eta4) {
338a7e14dcfSSatish Balay               /* Good step */
339a7e14dcfSSatish Balay               tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust);
34053506e15SBarry Smith             } else {
341a7e14dcfSSatish Balay               /* Very good step */
342a7e14dcfSSatish Balay               tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust);
343a7e14dcfSSatish Balay             }
344a7e14dcfSSatish Balay           }
345a7e14dcfSSatish Balay         }
346a7e14dcfSSatish Balay       }
34753506e15SBarry Smith     } else {
348a7e14dcfSSatish Balay       /* Get predicted reduction */
3499566063dSJacob Faibussowitsch       PetscCall(KSPCGGetObjFcn(tao->ksp,&prered));
350a7e14dcfSSatish Balay       if (prered >= 0.0) {
351a7e14dcfSSatish Balay         /* The predicted reduction has the wrong sign.  This cannot
352a7e14dcfSSatish Balay            happen in infinite precision arithmetic.  Step should
353a7e14dcfSSatish Balay            be rejected! */
354a7e14dcfSSatish Balay         tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
355a7e14dcfSSatish Balay         tr_reject = 1;
35653506e15SBarry Smith       } else {
3579566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->solution, tl->W));
3589566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tl->W, 1.0, tao->stepdirection));
3599566063dSJacob Faibussowitsch         PetscCall(TaoComputeObjective(tao, tl->W, &ftrial));
360a7e14dcfSSatish Balay         if (PetscIsInfOrNanReal(ftrial)) {
361a7e14dcfSSatish Balay           tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
362a7e14dcfSSatish Balay           tr_reject = 1;
36353506e15SBarry Smith         } else {
3649566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
365a7e14dcfSSatish Balay 
366a7e14dcfSSatish Balay           actred = f - ftrial;
367a7e14dcfSSatish Balay           prered = -prered;
368a7e14dcfSSatish Balay           if ((PetscAbsScalar(actred) <= tl->epsilon) &&
369a7e14dcfSSatish Balay               (PetscAbsScalar(prered) <= tl->epsilon)) {
370a7e14dcfSSatish Balay             kappa = 1.0;
37153506e15SBarry Smith           } else {
372a7e14dcfSSatish Balay             kappa = actred / prered;
373a7e14dcfSSatish Balay           }
374a7e14dcfSSatish Balay 
375a7e14dcfSSatish Balay           tau_1 = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred);
376a7e14dcfSSatish Balay           tau_2 = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred);
377a7e14dcfSSatish Balay           tau_min = PetscMin(tau_1, tau_2);
378a7e14dcfSSatish Balay           tau_max = PetscMax(tau_1, tau_2);
379a7e14dcfSSatish Balay 
380a7e14dcfSSatish Balay           if (kappa >= 1.0 - tl->mu1) {
381a7e14dcfSSatish Balay             /* Great agreement; accept step and update radius */
382a7e14dcfSSatish Balay             if (tau_max < 1.0) {
383a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d);
38453506e15SBarry Smith             } else if (tau_max > tl->gamma4) {
385a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d);
38653506e15SBarry Smith             } else {
387a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tau_max * norm_d);
388a7e14dcfSSatish Balay             }
38953506e15SBarry Smith           } else if (kappa >= 1.0 - tl->mu2) {
390a7e14dcfSSatish Balay             /* Good agreement */
391a7e14dcfSSatish Balay 
392a7e14dcfSSatish Balay             if (tau_max < tl->gamma2) {
393a7e14dcfSSatish Balay               tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d);
39453506e15SBarry Smith             } else if (tau_max > tl->gamma3) {
395a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d);
396a7e14dcfSSatish Balay             } else if (tau_max < 1.0) {
397a7e14dcfSSatish Balay               tao->trust = tau_max * PetscMin(tao->trust, norm_d);
39853506e15SBarry Smith             } else {
399a7e14dcfSSatish Balay               tao->trust = PetscMax(tao->trust, tau_max * norm_d);
400a7e14dcfSSatish Balay             }
40153506e15SBarry Smith           } else {
402a7e14dcfSSatish Balay             /* Not good agreement */
403a7e14dcfSSatish Balay             if (tau_min > 1.0) {
404a7e14dcfSSatish Balay               tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d);
40553506e15SBarry Smith             } else if (tau_max < tl->gamma1) {
406a7e14dcfSSatish Balay               tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
40753506e15SBarry Smith             } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) {
408a7e14dcfSSatish Balay               tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d);
40953506e15SBarry Smith             } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) {
410a7e14dcfSSatish Balay               tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
41153506e15SBarry Smith             } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) {
412a7e14dcfSSatish Balay               tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
41353506e15SBarry Smith             } else {
414a7e14dcfSSatish Balay               tao->trust = tau_max * PetscMin(tao->trust, norm_d);
415a7e14dcfSSatish Balay             }
416a7e14dcfSSatish Balay             tr_reject = 1;
417a7e14dcfSSatish Balay           }
418a7e14dcfSSatish Balay         }
419a7e14dcfSSatish Balay       }
420a7e14dcfSSatish Balay     }
421a7e14dcfSSatish Balay 
422a7e14dcfSSatish Balay     if (tr_reject) {
423a7e14dcfSSatish Balay       /* The trust-region constraints rejected the step.  Apply a linesearch.
424a7e14dcfSSatish Balay          Check for descent direction. */
4259566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
426a7e14dcfSSatish Balay       if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
427a7e14dcfSSatish Balay         /* Newton step is not descent or direction produced Inf or NaN */
428a7e14dcfSSatish Balay 
4290c51296cSAlp Dener         if (!tl->bfgs_pre) {
430a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and updated
431a7e14dcfSSatish Balay              Must use gradient direction in this case */
4329566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->gradient, tao->stepdirection));
4339566063dSJacob Faibussowitsch           PetscCall(VecScale(tao->stepdirection, -1.0));
434a7e14dcfSSatish Balay           ++tl->grad;
435a7e14dcfSSatish Balay           stepType = NTL_GRADIENT;
43653506e15SBarry Smith         } else {
437a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
4389566063dSJacob Faibussowitsch           PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
4399566063dSJacob Faibussowitsch           PetscCall(VecScale(tao->stepdirection, -1.0));
440a7e14dcfSSatish Balay 
441a7e14dcfSSatish Balay           /* Check for success (descent direction) */
4429566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
443a7e14dcfSSatish Balay           if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
444a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
445a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case because
446a7e14dcfSSatish Balay                the first solve produces the scaled gradient direction,
447a7e14dcfSSatish Balay                which is guaranteed to be descent */
448a7e14dcfSSatish Balay 
449a7e14dcfSSatish Balay             /* Use steepest descent direction (scaled) */
4509566063dSJacob Faibussowitsch             PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
4519566063dSJacob Faibussowitsch             PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
4529566063dSJacob Faibussowitsch             PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
4539566063dSJacob Faibussowitsch             PetscCall(VecScale(tao->stepdirection, -1.0));
454a7e14dcfSSatish Balay 
455a7e14dcfSSatish Balay             bfgsUpdates = 1;
4560c51296cSAlp Dener             ++tl->grad;
4570c51296cSAlp Dener             stepType = NTL_GRADIENT;
45853506e15SBarry Smith           } else {
4599566063dSJacob Faibussowitsch             PetscCall(MatLMVMGetUpdateCount(tl->M, &bfgsUpdates));
460a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
461a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
4620c51296cSAlp Dener               ++tl->grad;
4630c51296cSAlp Dener               stepType = NTL_GRADIENT;
46453506e15SBarry Smith             } else {
465a7e14dcfSSatish Balay               ++tl->bfgs;
466a7e14dcfSSatish Balay               stepType = NTL_BFGS;
467a7e14dcfSSatish Balay             }
468a7e14dcfSSatish Balay           }
469a7e14dcfSSatish Balay         }
47053506e15SBarry Smith       } else {
471a7e14dcfSSatish Balay         /* Computed Newton step is descent */
472a7e14dcfSSatish Balay         ++tl->newt;
473a7e14dcfSSatish Balay         stepType = NTL_NEWTON;
474a7e14dcfSSatish Balay       }
475a7e14dcfSSatish Balay 
476a7e14dcfSSatish Balay       /* Perform the linesearch */
477a7e14dcfSSatish Balay       fold = f;
4789566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->solution, tl->Xold));
4799566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->gradient, tl->Gold));
480a7e14dcfSSatish Balay 
481a7e14dcfSSatish Balay       step = 1.0;
4829566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason));
4839566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
484a7e14dcfSSatish Balay 
48553506e15SBarry Smith       while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) {      /* Linesearch failed */
486a7e14dcfSSatish Balay         /* Linesearch failed */
487a7e14dcfSSatish Balay         f = fold;
4889566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Xold, tao->solution));
4899566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Gold, tao->gradient));
490a7e14dcfSSatish Balay 
491a7e14dcfSSatish Balay         switch(stepType) {
492a7e14dcfSSatish Balay         case NTL_NEWTON:
493a7e14dcfSSatish Balay           /* Failed to obtain acceptable iterate with Newton step */
494a7e14dcfSSatish Balay 
4950c51296cSAlp Dener           if (tl->bfgs_pre) {
496a7e14dcfSSatish Balay             /* We don't have the bfgs matrix around and being updated
497a7e14dcfSSatish Balay                Must use gradient direction in this case */
4989566063dSJacob Faibussowitsch             PetscCall(VecCopy(tao->gradient, tao->stepdirection));
499a7e14dcfSSatish Balay             ++tl->grad;
500a7e14dcfSSatish Balay             stepType = NTL_GRADIENT;
50153506e15SBarry Smith           } else {
502a7e14dcfSSatish Balay             /* Attempt to use the BFGS direction */
5039566063dSJacob Faibussowitsch             PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
504a7e14dcfSSatish Balay 
505a7e14dcfSSatish Balay             /* Check for success (descent direction) */
5069566063dSJacob Faibussowitsch             PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
507a7e14dcfSSatish Balay             if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
508a7e14dcfSSatish Balay               /* BFGS direction is not descent or direction produced
509a7e14dcfSSatish Balay                  not a number.  We can assert bfgsUpdates > 1 in this case
510a7e14dcfSSatish Balay                  Use steepest descent direction (scaled) */
5119566063dSJacob Faibussowitsch               PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
5129566063dSJacob Faibussowitsch               PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
5139566063dSJacob Faibussowitsch               PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
514a7e14dcfSSatish Balay 
515a7e14dcfSSatish Balay               bfgsUpdates = 1;
5160c51296cSAlp Dener               ++tl->grad;
5170c51296cSAlp Dener               stepType = NTL_GRADIENT;
51853506e15SBarry Smith             } else {
5199566063dSJacob Faibussowitsch               PetscCall(MatLMVMGetUpdateCount(tl->M, &bfgsUpdates));
520a7e14dcfSSatish Balay               if (1 == bfgsUpdates) {
521a7e14dcfSSatish Balay                 /* The first BFGS direction is always the scaled gradient */
5220c51296cSAlp Dener                 ++tl->grad;
5230c51296cSAlp Dener                 stepType = NTL_GRADIENT;
52453506e15SBarry Smith               } else {
525a7e14dcfSSatish Balay                 ++tl->bfgs;
526a7e14dcfSSatish Balay                 stepType = NTL_BFGS;
527a7e14dcfSSatish Balay               }
528a7e14dcfSSatish Balay             }
529a7e14dcfSSatish Balay           }
530a7e14dcfSSatish Balay           break;
531a7e14dcfSSatish Balay 
532a7e14dcfSSatish Balay         case NTL_BFGS:
533a7e14dcfSSatish Balay           /* Can only enter if pc_type == NTL_PC_BFGS
534a7e14dcfSSatish Balay              Failed to obtain acceptable iterate with BFGS step
535a7e14dcfSSatish Balay              Attempt to use the scaled gradient direction */
5369566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(tl->M, PETSC_FALSE));
5379566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(tl->M, tao->solution, tao->gradient));
5389566063dSJacob Faibussowitsch           PetscCall(MatSolve(tl->M, tao->gradient, tao->stepdirection));
539a7e14dcfSSatish Balay 
540a7e14dcfSSatish Balay           bfgsUpdates = 1;
541a7e14dcfSSatish Balay           ++tl->grad;
542a7e14dcfSSatish Balay           stepType = NTL_GRADIENT;
543a7e14dcfSSatish Balay           break;
544a7e14dcfSSatish Balay         }
5459566063dSJacob Faibussowitsch         PetscCall(VecScale(tao->stepdirection, -1.0));
546a7e14dcfSSatish Balay 
547a7e14dcfSSatish Balay         /* This may be incorrect; linesearch has values for stepmax and stepmin
548a7e14dcfSSatish Balay            that should be reset. */
549a7e14dcfSSatish Balay         step = 1.0;
5509566063dSJacob Faibussowitsch         PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason));
5519566063dSJacob Faibussowitsch         PetscCall(TaoAddLineSearchCounts(tao));
552a7e14dcfSSatish Balay       }
553a7e14dcfSSatish Balay 
55453506e15SBarry Smith       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
555a7e14dcfSSatish Balay         /* Failed to find an improving point */
556a7e14dcfSSatish Balay         f = fold;
5579566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Xold, tao->solution));
5589566063dSJacob Faibussowitsch         PetscCall(VecCopy(tl->Gold, tao->gradient));
559a7e14dcfSSatish Balay         tao->trust = 0.0;
560a7e14dcfSSatish Balay         step = 0.0;
561a7e14dcfSSatish Balay         tao->reason = TAO_DIVERGED_LS_FAILURE;
562a7e14dcfSSatish Balay         break;
56353506e15SBarry Smith       } else if (stepType == NTL_NEWTON) {
564a7e14dcfSSatish Balay         if (step < tl->nu1) {
565a7e14dcfSSatish Balay           /* Very bad step taken; reduce radius */
566a7e14dcfSSatish Balay           tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust);
56753506e15SBarry Smith         } else if (step < tl->nu2) {
568a7e14dcfSSatish Balay           /* Reasonably bad step taken; reduce radius */
569a7e14dcfSSatish Balay           tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust);
57053506e15SBarry Smith         } else if (step < tl->nu3) {
571a7e14dcfSSatish Balay           /* Reasonable step was taken; leave radius alone */
572a7e14dcfSSatish Balay           if (tl->omega3 < 1.0) {
573a7e14dcfSSatish Balay             tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust);
57453506e15SBarry Smith           } else if (tl->omega3 > 1.0) {
575a7e14dcfSSatish Balay             tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust);
576a7e14dcfSSatish Balay           }
57753506e15SBarry Smith         } else if (step < tl->nu4) {
578a7e14dcfSSatish Balay           /* Full step taken; increase the radius */
579a7e14dcfSSatish Balay           tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust);
58053506e15SBarry Smith         } else {
581a7e14dcfSSatish Balay           /* More than full step taken; increase the radius */
582a7e14dcfSSatish Balay           tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust);
583a7e14dcfSSatish Balay         }
58453506e15SBarry Smith       } else {
585a7e14dcfSSatish Balay         /* Newton step was not good; reduce the radius */
586a7e14dcfSSatish Balay         tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust);
587a7e14dcfSSatish Balay       }
58853506e15SBarry Smith     } else {
589a7e14dcfSSatish Balay       /* Trust-region step is accepted */
5909566063dSJacob Faibussowitsch       PetscCall(VecCopy(tl->W, tao->solution));
591a7e14dcfSSatish Balay       f = ftrial;
5929566063dSJacob Faibussowitsch       PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
593a7e14dcfSSatish Balay       ++tl->ntrust;
594a7e14dcfSSatish Balay     }
595a7e14dcfSSatish Balay 
596a7e14dcfSSatish Balay     /* The radius may have been increased; modify if it is too large */
597a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, tl->max_radius);
598a7e14dcfSSatish Balay 
599e4cb33bbSBarry Smith     /* Check for converged */
6009566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
6013c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER,"User provided compute function generated Not-a-Number");
602a7e14dcfSSatish Balay     needH = 1;
603a7e14dcfSSatish Balay 
6049566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
6059566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
6069566063dSJacob Faibussowitsch     PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
607a7e14dcfSSatish Balay   }
608a7e14dcfSSatish Balay   PetscFunctionReturn(0);
609a7e14dcfSSatish Balay }
610a7e14dcfSSatish Balay 
611a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
612441846f8SBarry Smith static PetscErrorCode TaoSetUp_NTL(Tao tao)
613a7e14dcfSSatish Balay {
614a7e14dcfSSatish Balay   TAO_NTL        *tl = (TAO_NTL *)tao->data;
615a7e14dcfSSatish Balay 
616a7e14dcfSSatish Balay   PetscFunctionBegin;
6179566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
6189566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
6199566063dSJacob Faibussowitsch   if (!tl->W) PetscCall(VecDuplicate(tao->solution, &tl->W));
6209566063dSJacob Faibussowitsch   if (!tl->Xold) PetscCall(VecDuplicate(tao->solution, &tl->Xold));
6219566063dSJacob Faibussowitsch   if (!tl->Gold) PetscCall(VecDuplicate(tao->solution, &tl->Gold));
62283c8fe1dSLisandro Dalcin   tl->bfgs_pre = NULL;
62383c8fe1dSLisandro Dalcin   tl->M = NULL;
624a7e14dcfSSatish Balay   PetscFunctionReturn(0);
625a7e14dcfSSatish Balay }
626a7e14dcfSSatish Balay 
627a7e14dcfSSatish Balay /*------------------------------------------------------------*/
628441846f8SBarry Smith static PetscErrorCode TaoDestroy_NTL(Tao tao)
629a7e14dcfSSatish Balay {
630a7e14dcfSSatish Balay   TAO_NTL        *tl = (TAO_NTL *)tao->data;
631a7e14dcfSSatish Balay 
632a7e14dcfSSatish Balay   PetscFunctionBegin;
633a7e14dcfSSatish Balay   if (tao->setupcalled) {
6349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tl->W));
6359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tl->Xold));
6369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tl->Gold));
637a7e14dcfSSatish Balay   }
638*a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
6399566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
640a7e14dcfSSatish Balay   PetscFunctionReturn(0);
641a7e14dcfSSatish Balay }
642a7e14dcfSSatish Balay 
643a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6444416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_NTL(PetscOptionItems *PetscOptionsObject,Tao tao)
645a7e14dcfSSatish Balay {
646a7e14dcfSSatish Balay   TAO_NTL        *tl = (TAO_NTL *)tao->data;
647a7e14dcfSSatish Balay 
648a7e14dcfSSatish Balay   PetscFunctionBegin;
649d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Newton trust region with line search method for unconstrained optimization");
6509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_ntl_init_type", "radius initialization type", "", NTL_INIT, NTL_INIT_TYPES, NTL_INIT[tl->init_type], &tl->init_type,NULL));
6519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_ntl_update_type", "radius update type", "", NTL_UPDATE, NTL_UPDATE_TYPES, NTL_UPDATE[tl->update_type], &tl->update_type,NULL));
6529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1,NULL));
6539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2,NULL));
6549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3,NULL));
6559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4,NULL));
6569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1,NULL));
6579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2,NULL));
6589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3,NULL));
6599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4,NULL));
6609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5,NULL));
6619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1,NULL));
6629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2,NULL));
6639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3,NULL));
6649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4,NULL));
6659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1,NULL));
6669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2,NULL));
6679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3,NULL));
6689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4,NULL));
6699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5,NULL));
6709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i,NULL));
6719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i,NULL));
6729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i,NULL));
6739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i,NULL));
6749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i,NULL));
6759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i,NULL));
6769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i,NULL));
6779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1,NULL));
6789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2,NULL));
6799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1,NULL));
6809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2,NULL));
6819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3,NULL));
6829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4,NULL));
6839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta,NULL));
6849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius,NULL));
6859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius,NULL));
6869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon,NULL));
687d0609cedSBarry Smith   PetscOptionsHeadEnd();
6889566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
6899566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
690a7e14dcfSSatish Balay   PetscFunctionReturn(0);
691a7e14dcfSSatish Balay }
692a7e14dcfSSatish Balay 
693a7e14dcfSSatish Balay /*------------------------------------------------------------*/
694441846f8SBarry Smith static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer)
695a7e14dcfSSatish Balay {
696a7e14dcfSSatish Balay   TAO_NTL        *tl = (TAO_NTL *)tao->data;
697a7e14dcfSSatish Balay   PetscBool      isascii;
698a7e14dcfSSatish Balay 
699a7e14dcfSSatish Balay   PetscFunctionBegin;
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
701a7e14dcfSSatish Balay   if (isascii) {
7029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
70363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Trust-region steps: %" PetscInt_FMT "\n", tl->ntrust));
70463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton search steps: %" PetscInt_FMT "\n", tl->newt));
70563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS search steps: %" PetscInt_FMT "\n", tl->bfgs));
70663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient search steps: %" PetscInt_FMT "\n", tl->grad));
7079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
708a7e14dcfSSatish Balay   }
709a7e14dcfSSatish Balay   PetscFunctionReturn(0);
710a7e14dcfSSatish Balay }
711a7e14dcfSSatish Balay 
712a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
7131522df2eSJason Sarich /*MC
7143850be85SAlp Dener   TAONTL - Newton's method with trust region globalization and line search fallback.
7151522df2eSJason Sarich   At each iteration, the Newton trust region method solves the system for d
7161522df2eSJason Sarich   and performs a line search in the d direction:
7171522df2eSJason Sarich 
7181522df2eSJason Sarich             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
7191522df2eSJason Sarich 
7201522df2eSJason Sarich   Options Database Keys:
7219d0a60b2SAlp Dener + -tao_ntl_init_type - "constant","direction","interpolation"
7221522df2eSJason Sarich . -tao_ntl_update_type - "reduction","interpolation"
7231522df2eSJason Sarich . -tao_ntl_min_radius - lower bound on trust region radius
7241522df2eSJason Sarich . -tao_ntl_max_radius - upper bound on trust region radius
7251522df2eSJason Sarich . -tao_ntl_epsilon - tolerance for accepting actual / predicted reduction
7261522df2eSJason Sarich . -tao_ntl_mu1_i - mu1 interpolation init factor
7271522df2eSJason Sarich . -tao_ntl_mu2_i - mu2 interpolation init factor
7281522df2eSJason Sarich . -tao_ntl_gamma1_i - gamma1 interpolation init factor
7291522df2eSJason Sarich . -tao_ntl_gamma2_i - gamma2 interpolation init factor
7301522df2eSJason Sarich . -tao_ntl_gamma3_i - gamma3 interpolation init factor
7311522df2eSJason Sarich . -tao_ntl_gamma4_i - gamma4 interpolation init factor
7328966356dSPierre Jolivet . -tao_ntl_theta_i - theta1 interpolation init factor
7331522df2eSJason Sarich . -tao_ntl_eta1 - eta1 reduction update factor
7341522df2eSJason Sarich . -tao_ntl_eta2 - eta2 reduction update factor
7351522df2eSJason Sarich . -tao_ntl_eta3 - eta3 reduction update factor
7361522df2eSJason Sarich . -tao_ntl_eta4 - eta4 reduction update factor
7371522df2eSJason Sarich . -tao_ntl_alpha1 - alpha1 reduction update factor
7381522df2eSJason Sarich . -tao_ntl_alpha2 - alpha2 reduction update factor
7391522df2eSJason Sarich . -tao_ntl_alpha3 - alpha3 reduction update factor
7401522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor
7411522df2eSJason Sarich . -tao_ntl_alpha4 - alpha4 reduction update factor
7421522df2eSJason Sarich . -tao_ntl_mu1 - mu1 interpolation update
7431522df2eSJason Sarich . -tao_ntl_mu2 - mu2 interpolation update
7441522df2eSJason Sarich . -tao_ntl_gamma1 - gamma1 interpolcation update
7451522df2eSJason Sarich . -tao_ntl_gamma2 - gamma2 interpolcation update
7461522df2eSJason Sarich . -tao_ntl_gamma3 - gamma3 interpolcation update
7471522df2eSJason Sarich . -tao_ntl_gamma4 - gamma4 interpolation update
7481522df2eSJason Sarich - -tao_ntl_theta - theta1 interpolation update
7491522df2eSJason Sarich 
7501eb8069cSJason Sarich   Level: beginner
7511522df2eSJason Sarich M*/
752728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NTL(Tao tao)
753a7e14dcfSSatish Balay {
754a7e14dcfSSatish Balay   TAO_NTL        *tl;
7558caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
75653506e15SBarry Smith 
757a7e14dcfSSatish Balay   PetscFunctionBegin;
7589566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&tl));
759a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NTL;
760a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NTL;
761a7e14dcfSSatish Balay   tao->ops->view = TaoView_NTL;
762a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NTL;
763a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NTL;
764a7e14dcfSSatish Balay 
7656552cf8aSJason Sarich   /* Override default settings (unless already changed) */
7666552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
7676552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
7686552cf8aSJason Sarich 
769a7e14dcfSSatish Balay   tao->data = (void*)tl;
770a7e14dcfSSatish Balay 
771a7e14dcfSSatish Balay   /* Default values for trust-region radius update based on steplength */
772a7e14dcfSSatish Balay   tl->nu1 = 0.25;
773a7e14dcfSSatish Balay   tl->nu2 = 0.50;
774a7e14dcfSSatish Balay   tl->nu3 = 1.00;
775a7e14dcfSSatish Balay   tl->nu4 = 1.25;
776a7e14dcfSSatish Balay 
777a7e14dcfSSatish Balay   tl->omega1 = 0.25;
778a7e14dcfSSatish Balay   tl->omega2 = 0.50;
779a7e14dcfSSatish Balay   tl->omega3 = 1.00;
780a7e14dcfSSatish Balay   tl->omega4 = 2.00;
781a7e14dcfSSatish Balay   tl->omega5 = 4.00;
782a7e14dcfSSatish Balay 
783a7e14dcfSSatish Balay   /* Default values for trust-region radius update based on reduction */
784a7e14dcfSSatish Balay   tl->eta1 = 1.0e-4;
785a7e14dcfSSatish Balay   tl->eta2 = 0.25;
786a7e14dcfSSatish Balay   tl->eta3 = 0.50;
787a7e14dcfSSatish Balay   tl->eta4 = 0.90;
788a7e14dcfSSatish Balay 
789a7e14dcfSSatish Balay   tl->alpha1 = 0.25;
790a7e14dcfSSatish Balay   tl->alpha2 = 0.50;
791a7e14dcfSSatish Balay   tl->alpha3 = 1.00;
792a7e14dcfSSatish Balay   tl->alpha4 = 2.00;
793a7e14dcfSSatish Balay   tl->alpha5 = 4.00;
794a7e14dcfSSatish Balay 
795a7e14dcfSSatish Balay   /* Default values for trust-region radius update based on interpolation */
796a7e14dcfSSatish Balay   tl->mu1 = 0.10;
797a7e14dcfSSatish Balay   tl->mu2 = 0.50;
798a7e14dcfSSatish Balay 
799a7e14dcfSSatish Balay   tl->gamma1 = 0.25;
800a7e14dcfSSatish Balay   tl->gamma2 = 0.50;
801a7e14dcfSSatish Balay   tl->gamma3 = 2.00;
802a7e14dcfSSatish Balay   tl->gamma4 = 4.00;
803a7e14dcfSSatish Balay 
804a7e14dcfSSatish Balay   tl->theta = 0.05;
805a7e14dcfSSatish Balay 
806a7e14dcfSSatish Balay   /* Default values for trust region initialization based on interpolation */
807a7e14dcfSSatish Balay   tl->mu1_i = 0.35;
808a7e14dcfSSatish Balay   tl->mu2_i = 0.50;
809a7e14dcfSSatish Balay 
810a7e14dcfSSatish Balay   tl->gamma1_i = 0.0625;
811a7e14dcfSSatish Balay   tl->gamma2_i = 0.5;
812a7e14dcfSSatish Balay   tl->gamma3_i = 2.0;
813a7e14dcfSSatish Balay   tl->gamma4_i = 5.0;
814a7e14dcfSSatish Balay 
815a7e14dcfSSatish Balay   tl->theta_i = 0.25;
816a7e14dcfSSatish Balay 
817a7e14dcfSSatish Balay   /* Remaining parameters */
818a7e14dcfSSatish Balay   tl->min_radius = 1.0e-10;
819a7e14dcfSSatish Balay   tl->max_radius = 1.0e10;
820a7e14dcfSSatish Balay   tl->epsilon = 1.0e-6;
821a7e14dcfSSatish Balay 
822a7e14dcfSSatish Balay   tl->init_type       = NTL_INIT_INTERPOLATION;
823a7e14dcfSSatish Balay   tl->update_type     = NTL_UPDATE_REDUCTION;
824a7e14dcfSSatish Balay 
8259566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
8269566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch,(PetscObject)tao,1));
8279566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
8289566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
8299566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
8309566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp));
8319566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1));
8329566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix));
8339566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp,"tao_ntl_"));
8349566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp,KSPSTCG));
835a7e14dcfSSatish Balay   PetscFunctionReturn(0);
836a7e14dcfSSatish Balay }
837