xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
389371c9d4SSatish Balay static PetscErrorCode TaoSolve_NLS(Tao tao) {
39a7e14dcfSSatish Balay   TAO_NLS                     *nlsP = (TAO_NLS *)tao->data;
401daac38eSTodd Munson   KSPType                      ksp_type;
410ad3a497SAlp Dener   PetscBool                    is_nash, is_stcg, is_gltr, is_bfgs, is_jacobi, is_symmetric, sym_set;
42a7e14dcfSSatish Balay   KSPConvergedReason           ksp_reason;
431daac38eSTodd Munson   PC                           pc;
441daac38eSTodd Munson   TaoLineSearchConvergedReason ls_reason;
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   PetscReal fmin, ftrial, f_full, prered, actred, kappa, sigma;
47a7e14dcfSSatish Balay   PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
48a7e14dcfSSatish Balay   PetscReal f, fold, gdx, gnorm, pert;
49a7e14dcfSSatish Balay   PetscReal step   = 1.0;
50a7e14dcfSSatish Balay   PetscReal norm_d = 0.0, e_min;
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay   PetscInt stepType;
53a7e14dcfSSatish Balay   PetscInt bfgsUpdates = 0;
54a7e14dcfSSatish Balay   PetscInt n, N, kspits;
55b4690188SBarry Smith   PetscInt needH = 1;
56a7e14dcfSSatish Balay 
57a7e14dcfSSatish Balay   PetscInt i_max = 5;
58a7e14dcfSSatish Balay   PetscInt j_max = 1;
59a7e14dcfSSatish Balay   PetscInt i, j;
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay   PetscFunctionBegin;
62*48a46eb9SPierre 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"));
63a7e14dcfSSatish Balay 
64a7e14dcfSSatish Balay   /* Initialized variables */
65a7e14dcfSSatish Balay   pert = nlsP->sval;
66a7e14dcfSSatish Balay 
672d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
68a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
69a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
70a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
71a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
72a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
73a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
74a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
75a7e14dcfSSatish Balay 
76a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
77ba7fe8fbSTodd Munson      Command automatically ignored for other methods
78ba7fe8fbSTodd Munson      Will be reset during the first iteration
79ba7fe8fbSTodd Munson   */
809566063dSJacob Faibussowitsch   PetscCall(KSPGetType(tao->ksp, &ksp_type));
819566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPNASH, &is_nash));
829566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPSTCG, &is_stcg));
839566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type, KSPGLTR, &is_gltr));
841daac38eSTodd Munson 
859566063dSJacob Faibussowitsch   PetscCall(KSPCGSetRadius(tao->ksp, nlsP->max_radius));
86a7e14dcfSSatish Balay 
871daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
883c859ba3SBarry Smith     PetscCheck(tao->trust0 >= 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Initial radius negative");
89a7e14dcfSSatish Balay     tao->trust = tao->trust0;
90a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
91a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
92a7e14dcfSSatish Balay   }
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   /* Check convergence criteria */
959566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
969566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
973c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
98a7e14dcfSSatish Balay 
993ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1009566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
1019566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
102dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1033ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
104a7e14dcfSSatish Balay 
1050c51296cSAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
1069566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
1079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
1089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
1090c51296cSAlp Dener   if (is_bfgs) {
1100c51296cSAlp Dener     nlsP->bfgs_pre = pc;
1119566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(nlsP->bfgs_pre, &nlsP->M));
1129566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
1139566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
1149566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(nlsP->M, n, n, N, N));
1159566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(nlsP->M, tao->solution, tao->gradient));
1169566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(nlsP->M, &sym_set, &is_symmetric));
1173c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
1181baa6e33SBarry Smith   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
121a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
1221daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
123a7e14dcfSSatish Balay     switch (nlsP->init_type) {
124a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
125a7e14dcfSSatish Balay       /* Use the initial radius specified */
126a7e14dcfSSatish Balay       break;
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
129a7e14dcfSSatish Balay       /* Use the initial radius specified */
130a7e14dcfSSatish Balay       max_radius = 0.0;
131a7e14dcfSSatish Balay 
132a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
133a7e14dcfSSatish Balay         fmin  = f;
134a7e14dcfSSatish Balay         sigma = 0.0;
135a7e14dcfSSatish Balay 
136a7e14dcfSSatish Balay         if (needH) {
1379566063dSJacob Faibussowitsch           PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
138a7e14dcfSSatish Balay           needH = 0;
139a7e14dcfSSatish Balay         }
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
1429566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, nlsP->W));
1439566063dSJacob Faibussowitsch           PetscCall(VecAXPY(nlsP->W, -tao->trust / gnorm, tao->gradient));
1449566063dSJacob Faibussowitsch           PetscCall(TaoComputeObjective(tao, nlsP->W, &ftrial));
145a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
146a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
14787f595a5SBarry Smith           } else {
148a7e14dcfSSatish Balay             if (ftrial < fmin) {
149a7e14dcfSSatish Balay               fmin  = ftrial;
150a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
151a7e14dcfSSatish Balay             }
152a7e14dcfSSatish Balay 
1539566063dSJacob Faibussowitsch             PetscCall(MatMult(tao->hessian, tao->gradient, nlsP->D));
1549566063dSJacob Faibussowitsch             PetscCall(VecDot(tao->gradient, nlsP->D, &prered));
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
157a7e14dcfSSatish Balay             actred = f - ftrial;
15887f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
159a7e14dcfSSatish Balay               kappa = 1.0;
16087f595a5SBarry Smith             } else {
161a7e14dcfSSatish Balay               kappa = actred / prered;
162a7e14dcfSSatish Balay             }
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay             tau_1   = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
165a7e14dcfSSatish Balay             tau_2   = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
166a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
167a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
168a7e14dcfSSatish Balay 
16918cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu1_i) {
170a7e14dcfSSatish Balay               /* Great agreement */
171a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay               if (tau_max < 1.0) {
174a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
17587f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
176a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
17787f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
178a7e14dcfSSatish Balay                 tau = tau_1;
17987f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
180a7e14dcfSSatish Balay                 tau = tau_2;
18187f595a5SBarry Smith               } else {
182a7e14dcfSSatish Balay                 tau = tau_max;
183a7e14dcfSSatish Balay               }
18418cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu2_i) {
185a7e14dcfSSatish Balay               /* Good agreement */
186a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
187a7e14dcfSSatish Balay 
188a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
189a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
19087f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
191a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
19287f595a5SBarry Smith               } else {
193a7e14dcfSSatish Balay                 tau = tau_max;
194a7e14dcfSSatish Balay               }
19587f595a5SBarry Smith             } else {
196a7e14dcfSSatish Balay               /* Not good agreement */
197a7e14dcfSSatish Balay               if (tau_min > 1.0) {
198a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
19987f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
200a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
20187f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
202a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
20387f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
204a7e14dcfSSatish Balay                 tau = tau_1;
20587f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
206a7e14dcfSSatish Balay                 tau = tau_2;
20787f595a5SBarry Smith               } else {
208a7e14dcfSSatish Balay                 tau = tau_max;
209a7e14dcfSSatish Balay               }
210a7e14dcfSSatish Balay             }
211a7e14dcfSSatish Balay           }
212a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
213a7e14dcfSSatish Balay         }
214a7e14dcfSSatish Balay 
215a7e14dcfSSatish Balay         if (fmin < f) {
216a7e14dcfSSatish Balay           f = fmin;
2179566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
2189566063dSJacob Faibussowitsch           PetscCall(TaoComputeGradient(tao, tao->solution, tao->gradient));
219a7e14dcfSSatish Balay 
2209566063dSJacob Faibussowitsch           PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
2213c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute gradient generated Inf or NaN");
222a7e14dcfSSatish Balay           needH = 1;
223a7e14dcfSSatish Balay 
2249566063dSJacob Faibussowitsch           PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
2259566063dSJacob Faibussowitsch           PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
226dbbe0bcdSBarry Smith           PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
2273ecd9318SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
228a7e14dcfSSatish Balay         }
229a7e14dcfSSatish Balay       }
230a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
231a7e14dcfSSatish Balay 
232a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
233a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
234a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
235a7e14dcfSSatish Balay       break;
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay     default:
238a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
239a7e14dcfSSatish Balay       tao->trust = 0.0;
240a7e14dcfSSatish Balay       break;
241a7e14dcfSSatish Balay     }
242a7e14dcfSSatish Balay   }
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
245a7e14dcfSSatish Balay   nlsP->newt = 0;
246a7e14dcfSSatish Balay   nlsP->bfgs = 0;
247a7e14dcfSSatish Balay   nlsP->grad = 0;
248a7e14dcfSSatish Balay 
249a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2503ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
251e1e80dc8SAlp Dener     /* Call general purpose update function */
252dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
2538931d482SJason Sarich     ++tao->niter;
254ae93cb3cSJason Sarich     tao->ksp_its = 0;
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay     /* Compute the Hessian */
2571baa6e33SBarry Smith     if (needH) PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
258a7e14dcfSSatish Balay 
259a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
260a7e14dcfSSatish Balay     if (pert > 0) {
2619566063dSJacob Faibussowitsch       PetscCall(MatShift(tao->hessian, pert));
262*48a46eb9SPierre Jolivet       if (tao->hessian != tao->hessian_pre) PetscCall(MatShift(tao->hessian_pre, pert));
263a7e14dcfSSatish Balay     }
264a7e14dcfSSatish Balay 
2650c51296cSAlp Dener     if (nlsP->bfgs_pre) {
2669566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
267a7e14dcfSSatish Balay       ++bfgsUpdates;
268a7e14dcfSSatish Balay     }
269a7e14dcfSSatish Balay 
270a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
2719566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre));
2721daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
2739566063dSJacob Faibussowitsch       PetscCall(KSPCGSetRadius(tao->ksp, nlsP->max_radius));
2749566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
2759566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
276a7e14dcfSSatish Balay       tao->ksp_its += kspits;
277ae93cb3cSJason Sarich       tao->ksp_tot_its += kspits;
2789566063dSJacob Faibussowitsch       PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
279a7e14dcfSSatish Balay 
280a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
281a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
282a7e14dcfSSatish Balay         if (norm_d > 0.0) {
283a7e14dcfSSatish Balay           tao->trust = norm_d;
284a7e14dcfSSatish Balay 
285a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
286a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
287a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
28887f595a5SBarry Smith         } else {
289a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
290a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
291a7e14dcfSSatish Balay           tao->trust = tao->trust0;
292a7e14dcfSSatish Balay 
293a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
294a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
295a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
296a7e14dcfSSatish Balay 
2979566063dSJacob Faibussowitsch           PetscCall(KSPCGSetRadius(tao->ksp, nlsP->max_radius));
2989566063dSJacob Faibussowitsch           PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
2999566063dSJacob Faibussowitsch           PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
300a7e14dcfSSatish Balay           tao->ksp_its += kspits;
301ae93cb3cSJason Sarich           tao->ksp_tot_its += kspits;
3029566063dSJacob Faibussowitsch           PetscCall(KSPCGGetNormD(tao->ksp, &norm_d));
303a7e14dcfSSatish Balay 
3043c859ba3SBarry Smith           PetscCheck(norm_d != 0.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Initial direction zero");
305a7e14dcfSSatish Balay         }
306a7e14dcfSSatish Balay       }
30787f595a5SBarry Smith     } else {
3089566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
3099566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
310a7e14dcfSSatish Balay       tao->ksp_its += kspits;
311ae93cb3cSJason Sarich       tao->ksp_tot_its += kspits;
312a7e14dcfSSatish Balay     }
3139566063dSJacob Faibussowitsch     PetscCall(VecScale(nlsP->D, -1.0));
3149566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(tao->ksp, &ksp_reason));
3150c51296cSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (nlsP->bfgs_pre)) {
316a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
317a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
3189566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
3199566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
320a7e14dcfSSatish Balay       bfgsUpdates = 1;
321a7e14dcfSSatish Balay     }
322a7e14dcfSSatish Balay 
323a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
324a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
32587f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
326a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
32787f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
328a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
32987f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
330a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
33187f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
332a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
33387f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
334a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
33587f595a5SBarry Smith     } else {
336a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
337a7e14dcfSSatish Balay     }
338a7e14dcfSSatish Balay 
339a7e14dcfSSatish Balay     /* Check for success (descent direction) */
3409566063dSJacob Faibussowitsch     PetscCall(VecDot(nlsP->D, tao->gradient, &gdx));
341a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
342a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
343a7e14dcfSSatish Balay          Update the perturbation for next time */
344a7e14dcfSSatish Balay       if (pert <= 0.0) {
345a7e14dcfSSatish Balay         /* Initialize the perturbation */
346a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
3471daac38eSTodd Munson         if (is_gltr) {
3489566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
349a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
350a7e14dcfSSatish Balay         }
35187f595a5SBarry Smith       } else {
352a7e14dcfSSatish Balay         /* Increase the perturbation */
353a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
354a7e14dcfSSatish Balay       }
355a7e14dcfSSatish Balay 
3560c51296cSAlp Dener       if (!nlsP->bfgs_pre) {
357a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
358a7e14dcfSSatish Balay            Must use gradient direction in this case */
3599566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->gradient, nlsP->D));
3609566063dSJacob Faibussowitsch         PetscCall(VecScale(nlsP->D, -1.0));
361a7e14dcfSSatish Balay         ++nlsP->grad;
362a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
36387f595a5SBarry Smith       } else {
364a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
3659566063dSJacob Faibussowitsch         PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
3669566063dSJacob Faibussowitsch         PetscCall(VecScale(nlsP->D, -1.0));
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay         /* Check for success (descent direction) */
3699566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, nlsP->D, &gdx));
370a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
371a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
372a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
373a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
374a7e14dcfSSatish Balay              which is guaranteed to be descent */
375a7e14dcfSSatish Balay 
376a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
3779566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
3789566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
3799566063dSJacob Faibussowitsch           PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
3809566063dSJacob Faibussowitsch           PetscCall(VecScale(nlsP->D, -1.0));
381a7e14dcfSSatish Balay 
382a7e14dcfSSatish Balay           bfgsUpdates = 1;
3830c51296cSAlp Dener           ++nlsP->grad;
3840c51296cSAlp Dener           stepType = NLS_GRADIENT;
38587f595a5SBarry Smith         } else {
3869566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(nlsP->M, &bfgsUpdates));
387a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
388a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
3890c51296cSAlp Dener             ++nlsP->grad;
3900c51296cSAlp Dener             stepType = NLS_GRADIENT;
39187f595a5SBarry Smith           } else {
392a7e14dcfSSatish Balay             ++nlsP->bfgs;
393a7e14dcfSSatish Balay             stepType = NLS_BFGS;
394a7e14dcfSSatish Balay           }
395a7e14dcfSSatish Balay         }
396a7e14dcfSSatish Balay       }
39787f595a5SBarry Smith     } else {
398a7e14dcfSSatish Balay       /* Computed Newton step is descent */
399a7e14dcfSSatish Balay       switch (ksp_reason) {
400a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
401a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
402a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
403a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
404a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
405a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
406a7e14dcfSSatish Balay         if (pert <= 0.0) {
407a7e14dcfSSatish Balay           /* Initialize the perturbation */
408a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4091daac38eSTodd Munson           if (is_gltr) {
4109566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
411a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
412a7e14dcfSSatish Balay           }
41387f595a5SBarry Smith         } else {
414a7e14dcfSSatish Balay           /* Increase the perturbation */
415a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
416a7e14dcfSSatish Balay         }
417a7e14dcfSSatish Balay         break;
418a7e14dcfSSatish Balay 
419a7e14dcfSSatish Balay       default:
420a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
421a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
4229371c9d4SSatish Balay         if (pert < nlsP->pmin) { pert = 0.0; }
423a7e14dcfSSatish Balay         break;
424a7e14dcfSSatish Balay       }
425a7e14dcfSSatish Balay 
426a7e14dcfSSatish Balay       ++nlsP->newt;
427a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
428a7e14dcfSSatish Balay     }
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay     /* Perform the linesearch */
431a7e14dcfSSatish Balay     fold = f;
4329566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, nlsP->Xold));
4339566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, nlsP->Gold));
434a7e14dcfSSatish Balay 
4359566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason));
4369566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
437a7e14dcfSSatish Balay 
43887f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */
439a7e14dcfSSatish Balay       f = fold;
4409566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Xold, tao->solution));
4419566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Gold, tao->gradient));
442a7e14dcfSSatish Balay 
443a7e14dcfSSatish Balay       switch (stepType) {
444a7e14dcfSSatish Balay       case NLS_NEWTON:
445a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
446a7e14dcfSSatish Balay            Update the perturbation for next time */
447a7e14dcfSSatish Balay         if (pert <= 0.0) {
448a7e14dcfSSatish Balay           /* Initialize the perturbation */
449a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4501daac38eSTodd Munson           if (is_gltr) {
4519566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
452a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
453a7e14dcfSSatish Balay           }
45487f595a5SBarry Smith         } else {
455a7e14dcfSSatish Balay           /* Increase the perturbation */
456a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
457a7e14dcfSSatish Balay         }
458a7e14dcfSSatish Balay 
4590c51296cSAlp Dener         if (!nlsP->bfgs_pre) {
460a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
461a7e14dcfSSatish Balay              Must use gradient direction in this case */
4629566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->gradient, nlsP->D));
463a7e14dcfSSatish Balay           ++nlsP->grad;
464a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
46587f595a5SBarry Smith         } else {
466a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
4679566063dSJacob Faibussowitsch           PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
468a7e14dcfSSatish Balay           /* Check for success (descent direction) */
4699566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->solution, nlsP->D, &gdx));
470a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
471a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
472a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
473a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
4749566063dSJacob Faibussowitsch             PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
4759566063dSJacob Faibussowitsch             PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
4769566063dSJacob Faibussowitsch             PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
477a7e14dcfSSatish Balay 
478a7e14dcfSSatish Balay             bfgsUpdates = 1;
4790c51296cSAlp Dener             ++nlsP->grad;
4800c51296cSAlp Dener             stepType = NLS_GRADIENT;
48187f595a5SBarry Smith           } else {
482a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
483a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
4840c51296cSAlp Dener               ++nlsP->grad;
4850c51296cSAlp Dener               stepType = NLS_GRADIENT;
48687f595a5SBarry Smith             } else {
487a7e14dcfSSatish Balay               ++nlsP->bfgs;
488a7e14dcfSSatish Balay               stepType = NLS_BFGS;
489a7e14dcfSSatish Balay             }
490a7e14dcfSSatish Balay           }
491a7e14dcfSSatish Balay         }
492a7e14dcfSSatish Balay         break;
493a7e14dcfSSatish Balay 
494a7e14dcfSSatish Balay       case NLS_BFGS:
495a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
496a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
497a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
4989566063dSJacob Faibussowitsch         PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
4999566063dSJacob Faibussowitsch         PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
5009566063dSJacob Faibussowitsch         PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
501a7e14dcfSSatish Balay 
502a7e14dcfSSatish Balay         bfgsUpdates = 1;
503a7e14dcfSSatish Balay         ++nlsP->grad;
504a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
505a7e14dcfSSatish Balay         break;
506a7e14dcfSSatish Balay       }
5079566063dSJacob Faibussowitsch       PetscCall(VecScale(nlsP->D, -1.0));
508a7e14dcfSSatish Balay 
5099566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason));
5109566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full));
5119566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
512a7e14dcfSSatish Balay     }
513a7e14dcfSSatish Balay 
51487f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
515a7e14dcfSSatish Balay       /* Failed to find an improving point */
516a7e14dcfSSatish Balay       f = fold;
5179566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Xold, tao->solution));
5189566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Gold, tao->gradient));
519a7e14dcfSSatish Balay       step        = 0.0;
520a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
521a7e14dcfSSatish Balay       break;
522a7e14dcfSSatish Balay     }
523a7e14dcfSSatish Balay 
524a7e14dcfSSatish Balay     /* Update trust region radius */
5251daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
526a7e14dcfSSatish Balay       switch (nlsP->update_type) {
527a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
528a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
529a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
530a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
531a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
53287f595a5SBarry Smith           } else if (step < nlsP->nu2) {
533a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
534a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
53587f595a5SBarry Smith           } else if (step < nlsP->nu3) {
536a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
537a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
538a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
53987f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
540a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
541a7e14dcfSSatish Balay             }
54287f595a5SBarry Smith           } else if (step < nlsP->nu4) {
543a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
544a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
54587f595a5SBarry Smith           } else {
546a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
547a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
548a7e14dcfSSatish Balay           }
54987f595a5SBarry Smith         } else {
550a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
551a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
552a7e14dcfSSatish Balay         }
553a7e14dcfSSatish Balay         break;
554a7e14dcfSSatish Balay 
555a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
556a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
557a7e14dcfSSatish Balay           /*  Get predicted reduction */
558a7e14dcfSSatish Balay 
5599566063dSJacob Faibussowitsch           PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
560a7e14dcfSSatish Balay           if (prered >= 0.0) {
561a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
562a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
563a7e14dcfSSatish Balay             /*  be rejected! */
564a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
56587f595a5SBarry Smith           } else {
566a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
567a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
56887f595a5SBarry Smith             } else {
569a7e14dcfSSatish Balay               /*  Compute and actual reduction */
570a7e14dcfSSatish Balay               actred = fold - f_full;
571a7e14dcfSSatish Balay               prered = -prered;
5729371c9d4SSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
573a7e14dcfSSatish Balay                 kappa = 1.0;
57487f595a5SBarry Smith               } else {
575a7e14dcfSSatish Balay                 kappa = actred / prered;
576a7e14dcfSSatish Balay               }
577a7e14dcfSSatish Balay 
578a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
579a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
580a7e14dcfSSatish Balay                 /*  Very bad step */
581a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
58287f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
583a7e14dcfSSatish Balay                 /*  Marginal bad step */
584a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
58587f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
586a7e14dcfSSatish Balay                 /*  Reasonable step */
587a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
588a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
58987f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
590a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
591a7e14dcfSSatish Balay                 }
59287f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
593a7e14dcfSSatish Balay                 /*  Good step */
594a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
59587f595a5SBarry Smith               } else {
596a7e14dcfSSatish Balay                 /*  Very good step */
597a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
598a7e14dcfSSatish Balay               }
599a7e14dcfSSatish Balay             }
600a7e14dcfSSatish Balay           }
60187f595a5SBarry Smith         } else {
602a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
603a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
604a7e14dcfSSatish Balay         }
605a7e14dcfSSatish Balay         break;
606a7e14dcfSSatish Balay 
607a7e14dcfSSatish Balay       default:
608a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
6099566063dSJacob Faibussowitsch           PetscCall(KSPCGGetObjFcn(tao->ksp, &prered));
610a7e14dcfSSatish Balay           if (prered >= 0.0) {
611a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
612a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
613a7e14dcfSSatish Balay             /*  be rejected! */
614a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
61587f595a5SBarry Smith           } else {
616a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
617a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
61887f595a5SBarry Smith             } else {
619a7e14dcfSSatish Balay               actred = fold - f_full;
620a7e14dcfSSatish Balay               prered = -prered;
62187f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
622a7e14dcfSSatish Balay                 kappa = 1.0;
62387f595a5SBarry Smith               } else {
624a7e14dcfSSatish Balay                 kappa = actred / prered;
625a7e14dcfSSatish Balay               }
626a7e14dcfSSatish Balay 
627a7e14dcfSSatish Balay               tau_1   = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
628a7e14dcfSSatish Balay               tau_2   = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
629a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
630a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
631a7e14dcfSSatish Balay 
632a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
633a7e14dcfSSatish Balay                 /*  Great agreement */
634a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
635a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
63687f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
637a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
63887f595a5SBarry Smith                 } else {
639a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
640a7e14dcfSSatish Balay                 }
64187f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
642a7e14dcfSSatish Balay                 /*  Good agreement */
643a7e14dcfSSatish Balay 
644a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
645a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
64687f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
647a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
64887f595a5SBarry Smith                 } else if (tau_max < 1.0) {
649a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
65087f595a5SBarry Smith                 } else {
651a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
652a7e14dcfSSatish Balay                 }
65387f595a5SBarry Smith               } else {
654a7e14dcfSSatish Balay                 /*  Not good agreement */
655a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
656a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
65787f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
658a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
65987f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
660a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
66187f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
662a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
66387f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
664a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
66587f595a5SBarry Smith                 } else {
666a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
667a7e14dcfSSatish Balay                 }
668a7e14dcfSSatish Balay               }
669a7e14dcfSSatish Balay             }
670a7e14dcfSSatish Balay           }
67187f595a5SBarry Smith         } else {
672a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
673a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
674a7e14dcfSSatish Balay         }
675a7e14dcfSSatish Balay         break;
676a7e14dcfSSatish Balay       }
677a7e14dcfSSatish Balay 
678a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
679a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
680a7e14dcfSSatish Balay     }
681a7e14dcfSSatish Balay 
682a7e14dcfSSatish Balay     /*  Check for termination */
6839566063dSJacob Faibussowitsch     PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm));
6843c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Not-a-Number");
685a7e14dcfSSatish Balay     needH = 1;
6869566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
6879566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
688dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
689a7e14dcfSSatish Balay   }
690a7e14dcfSSatish Balay   PetscFunctionReturn(0);
691a7e14dcfSSatish Balay }
692a7e14dcfSSatish Balay 
693a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
6949371c9d4SSatish Balay static PetscErrorCode TaoSetUp_NLS(Tao tao) {
695a7e14dcfSSatish Balay   TAO_NLS *nlsP = (TAO_NLS *)tao->data;
696a7e14dcfSSatish Balay 
697a7e14dcfSSatish Balay   PetscFunctionBegin;
6989566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
6999566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
7009566063dSJacob Faibussowitsch   if (!nlsP->W) PetscCall(VecDuplicate(tao->solution, &nlsP->W));
7019566063dSJacob Faibussowitsch   if (!nlsP->D) PetscCall(VecDuplicate(tao->solution, &nlsP->D));
7029566063dSJacob Faibussowitsch   if (!nlsP->Xold) PetscCall(VecDuplicate(tao->solution, &nlsP->Xold));
7039566063dSJacob Faibussowitsch   if (!nlsP->Gold) PetscCall(VecDuplicate(tao->solution, &nlsP->Gold));
70483c8fe1dSLisandro Dalcin   nlsP->M        = NULL;
70583c8fe1dSLisandro Dalcin   nlsP->bfgs_pre = NULL;
706a7e14dcfSSatish Balay   PetscFunctionReturn(0);
707a7e14dcfSSatish Balay }
708a7e14dcfSSatish Balay 
709a7e14dcfSSatish Balay /*------------------------------------------------------------*/
7109371c9d4SSatish Balay static PetscErrorCode TaoDestroy_NLS(Tao tao) {
711a7e14dcfSSatish Balay   TAO_NLS *nlsP = (TAO_NLS *)tao->data;
712a7e14dcfSSatish Balay 
713a7e14dcfSSatish Balay   PetscFunctionBegin;
714a7e14dcfSSatish Balay   if (tao->setupcalled) {
7159566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->D));
7169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->W));
7179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->Xold));
7189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->Gold));
719a7e14dcfSSatish Balay   }
720a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
7219566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
722a7e14dcfSSatish Balay   PetscFunctionReturn(0);
723a7e14dcfSSatish Balay }
724a7e14dcfSSatish Balay 
725a7e14dcfSSatish Balay /*------------------------------------------------------------*/
7269371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_NLS(Tao tao, PetscOptionItems *PetscOptionsObject) {
727a7e14dcfSSatish Balay   TAO_NLS *nlsP = (TAO_NLS *)tao->data;
728a7e14dcfSSatish Balay 
729a7e14dcfSSatish Balay   PetscFunctionBegin;
730d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton line search method for unconstrained optimization");
7319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, NULL));
7329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, NULL));
7339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval, NULL));
7349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin, NULL));
7359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax, NULL));
7369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac, NULL));
7379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin, NULL));
7389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax, NULL));
7399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac, NULL));
7409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac, NULL));
7419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac, NULL));
7429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac, NULL));
7439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1, NULL));
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2, NULL));
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3, NULL));
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4, NULL));
7479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1, NULL));
7489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2, NULL));
7499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3, NULL));
7509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4, NULL));
7519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5, NULL));
7529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1, NULL));
7539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2, NULL));
7549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3, NULL));
7559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4, NULL));
7569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1, NULL));
7579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2, NULL));
7589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3, NULL));
7599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4, NULL));
7609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5, NULL));
7619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i, NULL));
7629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i, NULL));
7639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i, NULL));
7649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i, NULL));
7659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i, NULL));
7669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i, NULL));
7679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i, NULL));
7689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1, NULL));
7699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2, NULL));
7709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1, NULL));
7719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2, NULL));
7729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3, NULL));
7739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4, NULL));
7749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta, NULL));
7759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius, NULL));
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius, NULL));
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon, NULL));
778d0609cedSBarry Smith   PetscOptionsHeadEnd();
7799566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
7809566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
781a7e14dcfSSatish Balay   PetscFunctionReturn(0);
782a7e14dcfSSatish Balay }
783a7e14dcfSSatish Balay 
784a7e14dcfSSatish Balay /*------------------------------------------------------------*/
7859371c9d4SSatish Balay static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer) {
786a7e14dcfSSatish Balay   TAO_NLS  *nlsP = (TAO_NLS *)tao->data;
787a7e14dcfSSatish Balay   PetscBool isascii;
788a7e14dcfSSatish Balay 
789a7e14dcfSSatish Balay   PetscFunctionBegin;
7909566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
791a7e14dcfSSatish Balay   if (isascii) {
7929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
79363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", nlsP->newt));
79463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", nlsP->bfgs));
79563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", nlsP->grad));
796a7e14dcfSSatish Balay 
79763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp atol: %" PetscInt_FMT "\n", nlsP->ksp_atol));
79863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %" PetscInt_FMT "\n", nlsP->ksp_rtol));
79963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %" PetscInt_FMT "\n", nlsP->ksp_ctol));
80063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp negc: %" PetscInt_FMT "\n", nlsP->ksp_negc));
80163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %" PetscInt_FMT "\n", nlsP->ksp_dtol));
80263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp iter: %" PetscInt_FMT "\n", nlsP->ksp_iter));
80363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp othr: %" PetscInt_FMT "\n", nlsP->ksp_othr));
8049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
805a7e14dcfSSatish Balay   }
806a7e14dcfSSatish Balay   PetscFunctionReturn(0);
807a7e14dcfSSatish Balay }
808a7e14dcfSSatish Balay 
809a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
8104aa34175SJason Sarich /*MC
8114aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
8124aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
8134aa34175SJason Sarich   system of equations to obtain the step diretion dk:
8144aa34175SJason Sarich               Hk dk = -gk
8154aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
8164aa34175SJason Sarich   solve
8174aa34175SJason Sarich               min_t f(xk + t d_k)
8184aa34175SJason Sarich 
8194aa34175SJason Sarich     Options Database Keys:
8209d0a60b2SAlp Dener + -tao_nls_init_type - "constant","direction","interpolation"
8214aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
8224aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
8234aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
8244aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
8254aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
8264aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
8274aa34175SJason Sarich . -tao_nls_pgfac - growth factor
8284aa34175SJason Sarich . -tao_nls_psfac - shrink factor
8294aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
8304aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
8314aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
8324aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
8334aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
8344aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
8354aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
8364aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
8374aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
8384aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
8394aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
8404aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
8414aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
8424aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
8434aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
8444aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
8454aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
8464aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
8474aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
8484aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
8494aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
8504aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
8514aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
8524aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
8531522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
8541522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
8551522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
8561522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
8571522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
8581522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
8591522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
8601eb8069cSJason Sarich 
8611eb8069cSJason Sarich   Level: beginner
8621522df2eSJason Sarich M*/
8634aa34175SJason Sarich 
8649371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao) {
865a7e14dcfSSatish Balay   TAO_NLS    *nlsP;
8668caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
867a7e14dcfSSatish Balay 
868a7e14dcfSSatish Balay   PetscFunctionBegin;
8699566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao, &nlsP));
870a7e14dcfSSatish Balay 
871a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_NLS;
872a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_NLS;
873a7e14dcfSSatish Balay   tao->ops->view           = TaoView_NLS;
874a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
875a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_NLS;
876a7e14dcfSSatish Balay 
8776552cf8aSJason Sarich   /* Override default settings (unless already changed) */
8786552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
8796552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
8806552cf8aSJason Sarich 
881a7e14dcfSSatish Balay   tao->data = (void *)nlsP;
882a7e14dcfSSatish Balay 
883a7e14dcfSSatish Balay   nlsP->sval  = 0.0;
884a7e14dcfSSatish Balay   nlsP->imin  = 1.0e-4;
885a7e14dcfSSatish Balay   nlsP->imax  = 1.0e+2;
886a7e14dcfSSatish Balay   nlsP->imfac = 1.0e-1;
887a7e14dcfSSatish Balay 
888a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
889a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
890a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
891a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
892a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
893a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
894a7e14dcfSSatish Balay 
895a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
896a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
897a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
898a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
899a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
900a7e14dcfSSatish Balay 
901a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
902a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
903a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
904a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
905a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
906a7e14dcfSSatish Balay 
907a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
908a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
909a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
910a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
911a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
912a7e14dcfSSatish Balay 
913a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
914a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
915a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
916a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
917a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
918a7e14dcfSSatish Balay 
919a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
920a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
921a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
922a7e14dcfSSatish Balay 
923a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
924a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
925a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
926a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
927a7e14dcfSSatish Balay 
928a7e14dcfSSatish Balay   nlsP->theta = 0.05;
929a7e14dcfSSatish Balay 
930a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
931a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
932a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
933a7e14dcfSSatish Balay 
934a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
935a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
936a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
937a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
938a7e14dcfSSatish Balay 
939a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
940a7e14dcfSSatish Balay 
941a7e14dcfSSatish Balay   /*  Remaining parameters */
942a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
943a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
944a7e14dcfSSatish Balay   nlsP->epsilon    = 1.0e-6;
945a7e14dcfSSatish Balay 
946a7e14dcfSSatish Balay   nlsP->init_type   = NLS_INIT_INTERPOLATION;
947a7e14dcfSSatish Balay   nlsP->update_type = NLS_UPDATE_STEP;
948a7e14dcfSSatish Balay 
9499566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
9509566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
9519566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
9529566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
9539566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
954a7e14dcfSSatish Balay 
955a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
9569566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
9579566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
9589566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
9599566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_nls_"));
9609566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
961a7e14dcfSSatish Balay   PetscFunctionReturn(0);
962a7e14dcfSSatish Balay }
963