xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision 5a884c48ab0c46bab83cd9bb8710f380fa6d8bcf)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
21daac38eSTodd Munson #include <../src/tao/unconstrained/impls/nls/nlsimpl.h>
3a7e14dcfSSatish Balay 
4aaa7dc30SBarry Smith #include <petscksp.h>
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay #define NLS_INIT_CONSTANT      0
7a7e14dcfSSatish Balay #define NLS_INIT_DIRECTION     1
8a7e14dcfSSatish Balay #define NLS_INIT_INTERPOLATION 2
9a7e14dcfSSatish Balay #define NLS_INIT_TYPES         3
10a7e14dcfSSatish Balay 
11a7e14dcfSSatish Balay #define NLS_UPDATE_STEP          0
12a7e14dcfSSatish Balay #define NLS_UPDATE_REDUCTION     1
13a7e14dcfSSatish Balay #define NLS_UPDATE_INTERPOLATION 2
14a7e14dcfSSatish Balay #define NLS_UPDATE_TYPES         3
15a7e14dcfSSatish Balay 
1687f595a5SBarry Smith static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"};
17a7e14dcfSSatish Balay 
1887f595a5SBarry Smith static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"};
19a7e14dcfSSatish Balay 
201daac38eSTodd Munson /*
21a7e14dcfSSatish Balay  Implements Newton's Method with a line search approach for solving
22a7e14dcfSSatish Balay  unconstrained minimization problems.  A More'-Thuente line search
23a7e14dcfSSatish Balay  is used to guarantee that the bfgs preconditioner remains positive
24a7e14dcfSSatish Balay  definite.
25a7e14dcfSSatish Balay 
26a7e14dcfSSatish Balay  The method can shift the Hessian matrix.  The shifting procedure is
27a7e14dcfSSatish Balay  adapted from the PATH algorithm for solving complementarity
28a7e14dcfSSatish Balay  problems.
29a7e14dcfSSatish Balay 
30a7e14dcfSSatish Balay  The linear system solve should be done with a conjugate gradient
311daac38eSTodd Munson  method, although any method can be used.
321daac38eSTodd Munson */
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay #define NLS_NEWTON   0
35a7e14dcfSSatish Balay #define NLS_BFGS     1
360c51296cSAlp Dener #define NLS_GRADIENT 2
37a7e14dcfSSatish Balay 
TaoSolve_NLS(Tao tao)38d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_NLS(Tao tao)
39d71ae5a4SJacob Faibussowitsch {
40a7e14dcfSSatish Balay   TAO_NLS                     *nlsP = (TAO_NLS *)tao->data;
411daac38eSTodd Munson   KSPType                      ksp_type;
420ad3a497SAlp Dener   PetscBool                    is_nash, is_stcg, is_gltr, is_bfgs, is_jacobi, is_symmetric, sym_set;
43a7e14dcfSSatish Balay   KSPConvergedReason           ksp_reason;
441daac38eSTodd Munson   PC                           pc;
451daac38eSTodd Munson   TaoLineSearchConvergedReason ls_reason;
46a7e14dcfSSatish Balay 
47a7e14dcfSSatish Balay   PetscReal fmin, ftrial, f_full, prered, actred, kappa, sigma;
48a7e14dcfSSatish Balay   PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
49a7e14dcfSSatish Balay   PetscReal f, fold, gdx, gnorm, pert;
50a7e14dcfSSatish Balay   PetscReal step   = 1.0;
51a7e14dcfSSatish Balay   PetscReal norm_d = 0.0, e_min;
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay   PetscInt stepType;
54a7e14dcfSSatish Balay   PetscInt bfgsUpdates = 0;
55a7e14dcfSSatish Balay   PetscInt n, N, kspits;
56b4690188SBarry Smith   PetscInt needH = 1;
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay   PetscInt i_max = 5;
59a7e14dcfSSatish Balay   PetscInt j_max = 1;
60a7e14dcfSSatish Balay   PetscInt i, j;
61a7e14dcfSSatish Balay 
62a7e14dcfSSatish Balay   PetscFunctionBegin;
6348a46eb9SPierre Jolivet   if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by nls algorithm\n"));
64a7e14dcfSSatish Balay 
65a7e14dcfSSatish Balay   /* Initialized variables */
66a7e14dcfSSatish Balay   pert = nlsP->sval;
67a7e14dcfSSatish Balay 
682d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
69a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
70a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
71a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
72a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
73a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
74a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
75a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
78ba7fe8fbSTodd Munson      Command automatically ignored for other methods
79ba7fe8fbSTodd Munson      Will be reset during the first iteration
80ba7fe8fbSTodd Munson   */
819566063dSJacob Faibussowitsch   PetscCall(KSPGetType(tao->ksp, &ksp_type));
829566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPNASH, &is_nash));
839566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPSTCG, &is_stcg));
849566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPGLTR, &is_gltr));
851daac38eSTodd Munson 
869566063dSJacob Faibussowitsch   PetscCall(KSPCGSetRadius(tao->ksp, nlsP->max_radius));
87a7e14dcfSSatish Balay 
881daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
893c859ba3SBarry Smith     PetscCheck(tao->trust0 >= 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Initial radius negative");
90a7e14dcfSSatish Balay     tao->trust = tao->trust0;
91a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
92a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
93a7e14dcfSSatish Balay   }
94a7e14dcfSSatish Balay 
95a7e14dcfSSatish Balay   /* Check convergence criteria */
969566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
979566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
9876c63389SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
99a7e14dcfSSatish Balay 
1003ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1019566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
1029566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
103dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1043ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
105a7e14dcfSSatish Balay 
1060c51296cSAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
1079566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
1089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
1100c51296cSAlp Dener   if (is_bfgs) {
1110c51296cSAlp Dener     nlsP->bfgs_pre = pc;
1129566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(nlsP->bfgs_pre, &nlsP->M));
1139566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
1149566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
1159566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(nlsP->M, n, n, N, N));
1169566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(nlsP->M, tao->solution, tao->gradient));
1179566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(nlsP->M, &sym_set, &is_symmetric));
1183c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
1191baa6e33SBarry Smith   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
122a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
1231daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
124a7e14dcfSSatish Balay     switch (nlsP->init_type) {
125a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
126a7e14dcfSSatish Balay       /* Use the initial radius specified */
127a7e14dcfSSatish Balay       break;
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
130a7e14dcfSSatish Balay       /* Use the initial radius specified */
131a7e14dcfSSatish Balay       max_radius = 0.0;
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
134a7e14dcfSSatish Balay         fmin  = f;
135a7e14dcfSSatish Balay         sigma = 0.0;
136a7e14dcfSSatish Balay 
137a7e14dcfSSatish Balay         if (needH) {
1389566063dSJacob Faibussowitsch           PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
139a7e14dcfSSatish Balay           needH = 0;
140a7e14dcfSSatish Balay         }
141a7e14dcfSSatish Balay 
142a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
1439566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, nlsP->W));
1449566063dSJacob Faibussowitsch           PetscCall(VecAXPY(nlsP->W, -tao->trust / gnorm, tao->gradient));
1459566063dSJacob Faibussowitsch           PetscCall(TaoComputeObjective(tao, nlsP->W, &ftrial));
146a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
147a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
14887f595a5SBarry Smith           } else {
149a7e14dcfSSatish Balay             if (ftrial < fmin) {
150a7e14dcfSSatish Balay               fmin  = ftrial;
151a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
152a7e14dcfSSatish Balay             }
153a7e14dcfSSatish Balay 
1549566063dSJacob Faibussowitsch             PetscCall(MatMult(tao->hessian, tao->gradient, nlsP->D));
1559566063dSJacob Faibussowitsch             PetscCall(VecDot(tao->gradient, nlsP->D, &prered));
156a7e14dcfSSatish Balay 
157a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
158a7e14dcfSSatish Balay             actred = f - ftrial;
15987f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
160a7e14dcfSSatish Balay               kappa = 1.0;
16187f595a5SBarry Smith             } else {
162a7e14dcfSSatish Balay               kappa = actred / prered;
163a7e14dcfSSatish Balay             }
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay             tau_1   = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
166a7e14dcfSSatish Balay             tau_2   = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
167a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
168a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
169a7e14dcfSSatish Balay 
17018cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu1_i) {
171a7e14dcfSSatish Balay               /* Great agreement */
172a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay               if (tau_max < 1.0) {
175a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
17687f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
177a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
17887f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
179a7e14dcfSSatish Balay                 tau = tau_1;
18087f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
181a7e14dcfSSatish Balay                 tau = tau_2;
18287f595a5SBarry Smith               } else {
183a7e14dcfSSatish Balay                 tau = tau_max;
184a7e14dcfSSatish Balay               }
18518cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu2_i) {
186a7e14dcfSSatish Balay               /* Good agreement */
187a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
190a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
19187f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
192a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
19387f595a5SBarry Smith               } else {
194a7e14dcfSSatish Balay                 tau = tau_max;
195a7e14dcfSSatish Balay               }
19687f595a5SBarry Smith             } else {
197a7e14dcfSSatish Balay               /* Not good agreement */
198a7e14dcfSSatish Balay               if (tau_min > 1.0) {
199a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
20087f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
201a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
20287f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
203a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
20487f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
205a7e14dcfSSatish Balay                 tau = tau_1;
20687f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
207a7e14dcfSSatish Balay                 tau = tau_2;
20887f595a5SBarry Smith               } else {
209a7e14dcfSSatish Balay                 tau = tau_max;
210a7e14dcfSSatish Balay               }
211a7e14dcfSSatish Balay             }
212a7e14dcfSSatish Balay           }
213a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
214a7e14dcfSSatish Balay         }
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay         if (fmin < f) {
217a7e14dcfSSatish Balay           f = fmin;
2189566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
2199566063dSJacob Faibussowitsch           PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
220a7e14dcfSSatish Balay 
2219566063dSJacob Faibussowitsch           PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
22276c63389SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute gradient generated infinity or NaN");
223a7e14dcfSSatish Balay           needH = 1;
224a7e14dcfSSatish Balay 
2259566063dSJacob Faibussowitsch           PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
2269566063dSJacob Faibussowitsch           PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
227dbbe0bcdSBarry Smith           PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
2283ba16761SJacob Faibussowitsch           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
229a7e14dcfSSatish Balay         }
230a7e14dcfSSatish Balay       }
231a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
232a7e14dcfSSatish Balay 
233a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
234a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
235a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
236a7e14dcfSSatish Balay       break;
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay     default:
239a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
240a7e14dcfSSatish Balay       tao->trust = 0.0;
241a7e14dcfSSatish Balay       break;
242a7e14dcfSSatish Balay     }
243a7e14dcfSSatish Balay   }
244a7e14dcfSSatish Balay 
245a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
246a7e14dcfSSatish Balay   nlsP->newt = 0;
247a7e14dcfSSatish Balay   nlsP->bfgs = 0;
248a7e14dcfSSatish Balay   nlsP->grad = 0;
249a7e14dcfSSatish Balay 
250a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2513ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
252e1e80dc8SAlp Dener     /* Call general purpose update function */
253270bebe6SStefano Zampini     if (tao->ops->update) {
254270bebe6SStefano Zampini       PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
255270bebe6SStefano Zampini       PetscCall(TaoComputeObjective(tao, tao->solution, &f));
256270bebe6SStefano Zampini     }
2578931d482SJason Sarich     ++tao->niter;
258ae93cb3cSJason Sarich     tao->ksp_its = 0;
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay     /* Compute the Hessian */
2611baa6e33SBarry Smith     if (needH) PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
262a7e14dcfSSatish Balay 
263a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
264a7e14dcfSSatish Balay     if (pert > 0) {
2659566063dSJacob Faibussowitsch       PetscCall(MatShift(tao->hessian, pert));
26648a46eb9SPierre Jolivet       if (tao->hessian != tao->hessian_pre) PetscCall(MatShift(tao->hessian_pre, pert));
267a7e14dcfSSatish Balay     }
268a7e14dcfSSatish Balay 
2690c51296cSAlp Dener     if (nlsP->bfgs_pre) {
2709566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
271a7e14dcfSSatish Balay       ++bfgsUpdates;
272a7e14dcfSSatish Balay     }
273a7e14dcfSSatish Balay 
274a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
2759566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre));
2761daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
2779566063dSJacob Faibussowitsch       PetscCall(KSPCGSetRadius(tao->ksp, nlsP->max_radius));
2789566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
2799566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
280a7e14dcfSSatish Balay       tao->ksp_its += kspits;
281ae93cb3cSJason Sarich       tao->ksp_tot_its += kspits;
2829566063dSJacob Faibussowitsch       PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
283a7e14dcfSSatish Balay 
284a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
285a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
286a7e14dcfSSatish Balay         if (norm_d > 0.0) {
287a7e14dcfSSatish Balay           tao->trust = norm_d;
288a7e14dcfSSatish Balay 
289a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
290a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
291a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
29287f595a5SBarry Smith         } else {
293a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
294a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
295a7e14dcfSSatish Balay           tao->trust = tao->trust0;
296a7e14dcfSSatish Balay 
297a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
298a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
299a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
300a7e14dcfSSatish Balay 
3019566063dSJacob Faibussowitsch           PetscCall(KSPCGSetRadius(tao->ksp, nlsP->max_radius));
3029566063dSJacob Faibussowitsch           PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
3039566063dSJacob Faibussowitsch           PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
304a7e14dcfSSatish Balay           tao->ksp_its += kspits;
305ae93cb3cSJason Sarich           tao->ksp_tot_its += kspits;
3069566063dSJacob Faibussowitsch           PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
307a7e14dcfSSatish Balay 
3083c859ba3SBarry Smith           PetscCheck(norm_d != 0.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Initial direction zero");
309a7e14dcfSSatish Balay         }
310a7e14dcfSSatish Balay       }
31187f595a5SBarry Smith     } else {
3129566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
3139566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
314a7e14dcfSSatish Balay       tao->ksp_its += kspits;
315ae93cb3cSJason Sarich       tao->ksp_tot_its += kspits;
316a7e14dcfSSatish Balay     }
3179566063dSJacob Faibussowitsch     PetscCall(VecScale(nlsP->D, -1.0));
3189566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(tao->ksp, &ksp_reason));
3190c51296cSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (nlsP->bfgs_pre)) {
320a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
321a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
3229566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
3239566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
324a7e14dcfSSatish Balay       bfgsUpdates = 1;
325a7e14dcfSSatish Balay     }
326a7e14dcfSSatish Balay 
327a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
328a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
32987f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
330a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
3314a221d59SStefano Zampini     } else if (KSP_CONVERGED_STEP_LENGTH == ksp_reason) {
332a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
3334a221d59SStefano Zampini     } else if (KSP_CONVERGED_NEG_CURVE == ksp_reason) {
334a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
33587f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
336a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
33787f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
338a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
33987f595a5SBarry Smith     } else {
340a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
341a7e14dcfSSatish Balay     }
342a7e14dcfSSatish Balay 
343a7e14dcfSSatish Balay     /* Check for success (descent direction) */
3449566063dSJacob Faibussowitsch     PetscCall(VecDot(nlsP->D, tao->gradient, &gdx));
345a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
34676c63389SBarry Smith       /* Newton step is not descent or direction produced infinity or NaN
347a7e14dcfSSatish Balay          Update the perturbation for next time */
348a7e14dcfSSatish Balay       if (pert <= 0.0) {
349a7e14dcfSSatish Balay         /* Initialize the perturbation */
350a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
3511daac38eSTodd Munson         if (is_gltr) {
3529566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
353a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
354a7e14dcfSSatish Balay         }
35587f595a5SBarry Smith       } else {
356a7e14dcfSSatish Balay         /* Increase the perturbation */
357a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
358a7e14dcfSSatish Balay       }
359a7e14dcfSSatish Balay 
3600c51296cSAlp Dener       if (!nlsP->bfgs_pre) {
361a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
362a7e14dcfSSatish Balay            Must use gradient direction in this case */
3639566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->gradient, nlsP->D));
3649566063dSJacob Faibussowitsch         PetscCall(VecScale(nlsP->D, -1.0));
365a7e14dcfSSatish Balay         ++nlsP->grad;
366a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
36787f595a5SBarry Smith       } else {
368a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
3699566063dSJacob Faibussowitsch         PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
3709566063dSJacob Faibussowitsch         PetscCall(VecScale(nlsP->D, -1.0));
371a7e14dcfSSatish Balay 
372a7e14dcfSSatish Balay         /* Check for success (descent direction) */
3739566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, nlsP->D, &gdx));
374a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
375a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
376a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
377a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
378a7e14dcfSSatish Balay              which is guaranteed to be descent */
379a7e14dcfSSatish Balay 
380a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
3819566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
3829566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
3839566063dSJacob Faibussowitsch           PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
3849566063dSJacob Faibussowitsch           PetscCall(VecScale(nlsP->D, -1.0));
385a7e14dcfSSatish Balay 
386a7e14dcfSSatish Balay           bfgsUpdates = 1;
3870c51296cSAlp Dener           ++nlsP->grad;
3880c51296cSAlp Dener           stepType = NLS_GRADIENT;
38987f595a5SBarry Smith         } else {
3909566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(nlsP->M, &bfgsUpdates));
391a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
392a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
3930c51296cSAlp Dener             ++nlsP->grad;
3940c51296cSAlp Dener             stepType = NLS_GRADIENT;
39587f595a5SBarry Smith           } else {
396a7e14dcfSSatish Balay             ++nlsP->bfgs;
397a7e14dcfSSatish Balay             stepType = NLS_BFGS;
398a7e14dcfSSatish Balay           }
399a7e14dcfSSatish Balay         }
400a7e14dcfSSatish Balay       }
40187f595a5SBarry Smith     } else {
402a7e14dcfSSatish Balay       /* Computed Newton step is descent */
403a7e14dcfSSatish Balay       switch (ksp_reason) {
404a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
405a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
406a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
407a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
4084a221d59SStefano Zampini       case KSP_CONVERGED_NEG_CURVE:
409a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
410a7e14dcfSSatish Balay         if (pert <= 0.0) {
411a7e14dcfSSatish Balay           /* Initialize the perturbation */
412a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4131daac38eSTodd Munson           if (is_gltr) {
4149566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
415a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
416a7e14dcfSSatish Balay           }
41787f595a5SBarry Smith         } else {
418a7e14dcfSSatish Balay           /* Increase the perturbation */
419a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
420a7e14dcfSSatish Balay         }
421a7e14dcfSSatish Balay         break;
422a7e14dcfSSatish Balay 
423a7e14dcfSSatish Balay       default:
424a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
425a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
426ad540459SPierre Jolivet         if (pert < nlsP->pmin) pert = 0.0;
427a7e14dcfSSatish Balay         break;
428a7e14dcfSSatish Balay       }
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay       ++nlsP->newt;
431a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
432a7e14dcfSSatish Balay     }
433a7e14dcfSSatish Balay 
434a7e14dcfSSatish Balay     /* Perform the linesearch */
435a7e14dcfSSatish Balay     fold = f;
4369566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, nlsP->Xold));
4379566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, nlsP->Gold));
438a7e14dcfSSatish Balay 
4399566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason));
4409566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
441a7e14dcfSSatish Balay 
44287f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */
443a7e14dcfSSatish Balay       f = fold;
4449566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Xold, tao->solution));
4459566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Gold, tao->gradient));
446a7e14dcfSSatish Balay 
447a7e14dcfSSatish Balay       switch (stepType) {
448a7e14dcfSSatish Balay       case NLS_NEWTON:
449a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
450a7e14dcfSSatish Balay            Update the perturbation for next time */
451a7e14dcfSSatish Balay         if (pert <= 0.0) {
452a7e14dcfSSatish Balay           /* Initialize the perturbation */
453a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4541daac38eSTodd Munson           if (is_gltr) {
4559566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
456a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
457a7e14dcfSSatish Balay           }
45887f595a5SBarry Smith         } else {
459a7e14dcfSSatish Balay           /* Increase the perturbation */
460a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
461a7e14dcfSSatish Balay         }
462a7e14dcfSSatish Balay 
4630c51296cSAlp Dener         if (!nlsP->bfgs_pre) {
464a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
465a7e14dcfSSatish Balay              Must use gradient direction in this case */
4669566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->gradient, nlsP->D));
467a7e14dcfSSatish Balay           ++nlsP->grad;
468a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
46987f595a5SBarry Smith         } else {
470a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
4719566063dSJacob Faibussowitsch           PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
472a7e14dcfSSatish Balay           /* Check for success (descent direction) */
473*79813928SHan Liu           PetscCall(VecDot(tao->gradient, nlsP->D, &gdx));
474a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
475a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
476a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
477a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
4789566063dSJacob Faibussowitsch             PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
4799566063dSJacob Faibussowitsch             PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
4809566063dSJacob Faibussowitsch             PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
481a7e14dcfSSatish Balay 
482a7e14dcfSSatish Balay             bfgsUpdates = 1;
4830c51296cSAlp Dener             ++nlsP->grad;
4840c51296cSAlp Dener             stepType = NLS_GRADIENT;
48587f595a5SBarry Smith           } else {
486a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
487a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
4880c51296cSAlp Dener               ++nlsP->grad;
4890c51296cSAlp Dener               stepType = NLS_GRADIENT;
49087f595a5SBarry Smith             } else {
491a7e14dcfSSatish Balay               ++nlsP->bfgs;
492a7e14dcfSSatish Balay               stepType = NLS_BFGS;
493a7e14dcfSSatish Balay             }
494a7e14dcfSSatish Balay           }
495a7e14dcfSSatish Balay         }
496a7e14dcfSSatish Balay         break;
497a7e14dcfSSatish Balay 
498a7e14dcfSSatish Balay       case NLS_BFGS:
499a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
500a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
501a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
5029566063dSJacob Faibussowitsch         PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
5039566063dSJacob Faibussowitsch         PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
5049566063dSJacob Faibussowitsch         PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
505a7e14dcfSSatish Balay 
506a7e14dcfSSatish Balay         bfgsUpdates = 1;
507a7e14dcfSSatish Balay         ++nlsP->grad;
508a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
509a7e14dcfSSatish Balay         break;
510a7e14dcfSSatish Balay       }
5119566063dSJacob Faibussowitsch       PetscCall(VecScale(nlsP->D, -1.0));
512a7e14dcfSSatish Balay 
5139566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason));
5149566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full));
5159566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
516a7e14dcfSSatish Balay     }
517a7e14dcfSSatish Balay 
51887f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
519a7e14dcfSSatish Balay       /* Failed to find an improving point */
520a7e14dcfSSatish Balay       f = fold;
5219566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Xold, tao->solution));
5229566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Gold, tao->gradient));
523a7e14dcfSSatish Balay       step        = 0.0;
524a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
525a7e14dcfSSatish Balay       break;
526a7e14dcfSSatish Balay     }
527a7e14dcfSSatish Balay 
528a7e14dcfSSatish Balay     /* Update trust region radius */
5291daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
530a7e14dcfSSatish Balay       switch (nlsP->update_type) {
531a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
532a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
533a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
534a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
535a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
53687f595a5SBarry Smith           } else if (step < nlsP->nu2) {
537a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
538a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
53987f595a5SBarry Smith           } else if (step < nlsP->nu3) {
540a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
541a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
542a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
54387f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
544a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
545a7e14dcfSSatish Balay             }
54687f595a5SBarry Smith           } else if (step < nlsP->nu4) {
547a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
548a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
54987f595a5SBarry Smith           } else {
550a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
551a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
552a7e14dcfSSatish Balay           }
55387f595a5SBarry Smith         } else {
554a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
555a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
556a7e14dcfSSatish Balay         }
557a7e14dcfSSatish Balay         break;
558a7e14dcfSSatish Balay 
559a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
560a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
561a7e14dcfSSatish Balay           /*  Get predicted reduction */
562a7e14dcfSSatish Balay 
5639566063dSJacob Faibussowitsch           PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
564a7e14dcfSSatish Balay           if (prered >= 0.0) {
565a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
566a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
567a7e14dcfSSatish Balay             /*  be rejected! */
568a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
56987f595a5SBarry Smith           } else {
570a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
571a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
57287f595a5SBarry Smith             } else {
573a7e14dcfSSatish Balay               /*  Compute and actual reduction */
574a7e14dcfSSatish Balay               actred = fold - f_full;
575a7e14dcfSSatish Balay               prered = -prered;
5769371c9d4SSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
577a7e14dcfSSatish Balay                 kappa = 1.0;
57887f595a5SBarry Smith               } else {
579a7e14dcfSSatish Balay                 kappa = actred / prered;
580a7e14dcfSSatish Balay               }
581a7e14dcfSSatish Balay 
582a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
583a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
584a7e14dcfSSatish Balay                 /*  Very bad step */
585a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
58687f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
587a7e14dcfSSatish Balay                 /*  Marginal bad step */
588a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
58987f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
590a7e14dcfSSatish Balay                 /*  Reasonable step */
591a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
592a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
59387f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
594a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
595a7e14dcfSSatish Balay                 }
59687f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
597a7e14dcfSSatish Balay                 /*  Good step */
598a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
59987f595a5SBarry Smith               } else {
600a7e14dcfSSatish Balay                 /*  Very good step */
601a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
602a7e14dcfSSatish Balay               }
603a7e14dcfSSatish Balay             }
604a7e14dcfSSatish Balay           }
60587f595a5SBarry Smith         } else {
606a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
607a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
608a7e14dcfSSatish Balay         }
609a7e14dcfSSatish Balay         break;
610a7e14dcfSSatish Balay 
611a7e14dcfSSatish Balay       default:
612a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
6139566063dSJacob Faibussowitsch           PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
614a7e14dcfSSatish Balay           if (prered >= 0.0) {
615a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
616a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
617a7e14dcfSSatish Balay             /*  be rejected! */
618a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
61987f595a5SBarry Smith           } else {
620a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
621a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
62287f595a5SBarry Smith             } else {
623a7e14dcfSSatish Balay               actred = fold - f_full;
624a7e14dcfSSatish Balay               prered = -prered;
62587f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
626a7e14dcfSSatish Balay                 kappa = 1.0;
62787f595a5SBarry Smith               } else {
628a7e14dcfSSatish Balay                 kappa = actred / prered;
629a7e14dcfSSatish Balay               }
630a7e14dcfSSatish Balay 
631a7e14dcfSSatish Balay               tau_1   = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
632a7e14dcfSSatish Balay               tau_2   = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
633a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
634a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
635a7e14dcfSSatish Balay 
636a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
637a7e14dcfSSatish Balay                 /*  Great agreement */
638a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
639a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
64087f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
641a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
64287f595a5SBarry Smith                 } else {
643a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
644a7e14dcfSSatish Balay                 }
64587f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
646a7e14dcfSSatish Balay                 /*  Good agreement */
647a7e14dcfSSatish Balay 
648a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
649a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
65087f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
651a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
65287f595a5SBarry Smith                 } else if (tau_max < 1.0) {
653a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
65487f595a5SBarry Smith                 } else {
655a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
656a7e14dcfSSatish Balay                 }
65787f595a5SBarry Smith               } else {
658a7e14dcfSSatish Balay                 /*  Not good agreement */
659a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
660a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
66187f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
662a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
66387f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
664a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
66587f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
666a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
66787f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
668a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
66987f595a5SBarry Smith                 } else {
670a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
671a7e14dcfSSatish Balay                 }
672a7e14dcfSSatish Balay               }
673a7e14dcfSSatish Balay             }
674a7e14dcfSSatish Balay           }
67587f595a5SBarry Smith         } else {
676a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
677a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
678a7e14dcfSSatish Balay         }
679a7e14dcfSSatish Balay         break;
680a7e14dcfSSatish Balay       }
681a7e14dcfSSatish Balay 
682a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
683a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
684a7e14dcfSSatish Balay     }
685a7e14dcfSSatish Balay 
686a7e14dcfSSatish Balay     /*  Check for termination */
6879566063dSJacob Faibussowitsch     PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
6883c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Not-a-Number");
689a7e14dcfSSatish Balay     needH = 1;
6909566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
6919566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
692dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
693a7e14dcfSSatish Balay   }
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
695a7e14dcfSSatish Balay }
696a7e14dcfSSatish Balay 
TaoSetUp_NLS(Tao tao)697d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_NLS(Tao tao)
698d71ae5a4SJacob Faibussowitsch {
699a7e14dcfSSatish Balay   TAO_NLS *nlsP = (TAO_NLS *)tao->data;
700a7e14dcfSSatish Balay 
701a7e14dcfSSatish Balay   PetscFunctionBegin;
7029566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
7039566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
7049566063dSJacob Faibussowitsch   if (!nlsP->W) PetscCall(VecDuplicate(tao->solution, &nlsP->W));
7059566063dSJacob Faibussowitsch   if (!nlsP->D) PetscCall(VecDuplicate(tao->solution, &nlsP->D));
7069566063dSJacob Faibussowitsch   if (!nlsP->Xold) PetscCall(VecDuplicate(tao->solution, &nlsP->Xold));
7079566063dSJacob Faibussowitsch   if (!nlsP->Gold) PetscCall(VecDuplicate(tao->solution, &nlsP->Gold));
70883c8fe1dSLisandro Dalcin   nlsP->M        = NULL;
70983c8fe1dSLisandro Dalcin   nlsP->bfgs_pre = NULL;
7103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
711a7e14dcfSSatish Balay }
712a7e14dcfSSatish Balay 
TaoDestroy_NLS(Tao tao)713d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_NLS(Tao tao)
714d71ae5a4SJacob Faibussowitsch {
715a7e14dcfSSatish Balay   TAO_NLS *nlsP = (TAO_NLS *)tao->data;
716a7e14dcfSSatish Balay 
717a7e14dcfSSatish Balay   PetscFunctionBegin;
718a7e14dcfSSatish Balay   if (tao->setupcalled) {
7199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->D));
7209566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->W));
7219566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->Xold));
7229566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->Gold));
723a7e14dcfSSatish Balay   }
724a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
7259566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
727a7e14dcfSSatish Balay }
728a7e14dcfSSatish Balay 
TaoSetFromOptions_NLS(Tao tao,PetscOptionItems PetscOptionsObject)729ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(Tao tao, PetscOptionItems PetscOptionsObject)
730d71ae5a4SJacob Faibussowitsch {
731a7e14dcfSSatish Balay   TAO_NLS *nlsP = (TAO_NLS *)tao->data;
732a7e14dcfSSatish Balay 
733a7e14dcfSSatish Balay   PetscFunctionBegin;
734d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton line search method for unconstrained optimization");
7359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, NULL));
7369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, NULL));
7379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval, NULL));
7389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin, NULL));
7399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax, NULL));
7409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac, NULL));
7419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin, NULL));
7429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax, NULL));
7439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac, NULL));
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac, NULL));
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac, NULL));
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac, NULL));
7479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1, NULL));
7489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2, NULL));
7499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3, NULL));
7509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4, NULL));
7519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1, NULL));
7529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2, NULL));
7539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3, NULL));
7549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4, NULL));
7559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5, NULL));
7569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1, NULL));
7579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2, NULL));
7589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3, NULL));
7599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4, NULL));
7609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1, NULL));
7619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2, NULL));
7629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3, NULL));
7639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4, NULL));
7649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5, NULL));
7659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i, NULL));
7669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i, NULL));
7679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i, NULL));
7689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i, NULL));
7699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i, NULL));
7709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i, NULL));
7719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i, NULL));
7729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1, NULL));
7739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2, NULL));
7749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1, NULL));
7759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2, NULL));
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3, NULL));
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4, NULL));
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta, NULL));
7799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius, NULL));
7809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius, NULL));
7819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon, NULL));
782d0609cedSBarry Smith   PetscOptionsHeadEnd();
7839566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
7849566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
7853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
786a7e14dcfSSatish Balay }
787a7e14dcfSSatish Balay 
TaoView_NLS(Tao tao,PetscViewer viewer)788d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
789d71ae5a4SJacob Faibussowitsch {
790a7e14dcfSSatish Balay   TAO_NLS  *nlsP = (TAO_NLS *)tao->data;
791a7e14dcfSSatish Balay   PetscBool isascii;
792a7e14dcfSSatish Balay 
793a7e14dcfSSatish Balay   PetscFunctionBegin;
7949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
795a7e14dcfSSatish Balay   if (isascii) {
7969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
79763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", nlsP->newt));
79863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", nlsP->bfgs));
79963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", nlsP->grad));
800a7e14dcfSSatish Balay 
80163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp atol: %" PetscInt_FMT "\n", nlsP->ksp_atol));
80263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %" PetscInt_FMT "\n", nlsP->ksp_rtol));
80363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %" PetscInt_FMT "\n", nlsP->ksp_ctol));
80463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp negc: %" PetscInt_FMT "\n", nlsP->ksp_negc));
80563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %" PetscInt_FMT "\n", nlsP->ksp_dtol));
80663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp iter: %" PetscInt_FMT "\n", nlsP->ksp_iter));
80763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp othr: %" PetscInt_FMT "\n", nlsP->ksp_othr));
8089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
809a7e14dcfSSatish Balay   }
8103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
811a7e14dcfSSatish Balay }
812a7e14dcfSSatish Balay 
8134aa34175SJason Sarich /*MC
8144aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
8154aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
8161cb3f120SPierre Jolivet   system of equations to obtain the step direction dk:
8174aa34175SJason Sarich               Hk dk = -gk
8184aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
8194aa34175SJason Sarich   solve
8204aa34175SJason Sarich               min_t f(xk + t d_k)
8214aa34175SJason Sarich 
8224aa34175SJason Sarich     Options Database Keys:
8239d0a60b2SAlp Dener + -tao_nls_init_type - "constant","direction","interpolation"
8244aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
8254aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
8264aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
8274aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
8284aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
8294aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
8304aa34175SJason Sarich . -tao_nls_pgfac - growth factor
8314aa34175SJason Sarich . -tao_nls_psfac - shrink factor
8324aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
8334aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
8344aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
8354aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
8364aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
837da81f932SPierre Jolivet . -tao_nls_eta3 - good steplength; increase radius
8384aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
8394aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
8404aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
8414aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
8424aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
8434aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
8444aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
8454aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
8464aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
8474aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
8484aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
8494aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
8504aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
8514aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
8524aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
8534aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
8544aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
8554aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
8561522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
8571522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
8581522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
8591522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
8601522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
8611522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
8621522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
8631eb8069cSJason Sarich 
8641eb8069cSJason Sarich   Level: beginner
8651522df2eSJason Sarich M*/
8664aa34175SJason Sarich 
TaoCreate_NLS(Tao tao)867d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
868d71ae5a4SJacob Faibussowitsch {
869a7e14dcfSSatish Balay   TAO_NLS    *nlsP;
8708caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
871a7e14dcfSSatish Balay 
872a7e14dcfSSatish Balay   PetscFunctionBegin;
8734dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&nlsP));
874a7e14dcfSSatish Balay 
875a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_NLS;
876a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_NLS;
877a7e14dcfSSatish Balay   tao->ops->view           = TaoView_NLS;
878a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
879a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_NLS;
880a7e14dcfSSatish Balay 
8816552cf8aSJason Sarich   /* Override default settings (unless already changed) */
882606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
883606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 50);
884606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, trust0, 100.0);
8856552cf8aSJason Sarich 
886a7e14dcfSSatish Balay   tao->data = (void *)nlsP;
887a7e14dcfSSatish Balay 
888a7e14dcfSSatish Balay   nlsP->sval  = 0.0;
889a7e14dcfSSatish Balay   nlsP->imin  = 1.0e-4;
890a7e14dcfSSatish Balay   nlsP->imax  = 1.0e+2;
891a7e14dcfSSatish Balay   nlsP->imfac = 1.0e-1;
892a7e14dcfSSatish Balay 
893a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
894a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
895a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
896a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
897a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
898a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
899a7e14dcfSSatish Balay 
900a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
901a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
902a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
903a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
904a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
905a7e14dcfSSatish Balay 
906a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
907a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
908a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
909a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
910a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
911a7e14dcfSSatish Balay 
912a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
913a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
914a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
915a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
916a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
917a7e14dcfSSatish Balay 
918a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
919a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
920a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
921a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
922a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
923a7e14dcfSSatish Balay 
924a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
925a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
926a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
927a7e14dcfSSatish Balay 
928a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
929a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
930a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
931a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
932a7e14dcfSSatish Balay 
933a7e14dcfSSatish Balay   nlsP->theta = 0.05;
934a7e14dcfSSatish Balay 
935a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
936a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
937a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
938a7e14dcfSSatish Balay 
939a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
940a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
941a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
942a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
943a7e14dcfSSatish Balay 
944a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
945a7e14dcfSSatish Balay 
946a7e14dcfSSatish Balay   /*  Remaining parameters */
947a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
948a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
949a7e14dcfSSatish Balay   nlsP->epsilon    = 1.0e-6;
950a7e14dcfSSatish Balay 
951a7e14dcfSSatish Balay   nlsP->init_type   = NLS_INIT_INTERPOLATION;
952a7e14dcfSSatish Balay   nlsP->update_type = NLS_UPDATE_STEP;
953a7e14dcfSSatish Balay 
9549566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
9559566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
9569566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
9579566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
9589566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
959a7e14dcfSSatish Balay 
960a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
9619566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
9629566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
9639566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
9649566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_nls_"));
9659566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
9663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
967a7e14dcfSSatish Balay }
968