xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
38441846f8SBarry Smith static PetscErrorCode TaoSolve_NLS(Tao tao)
39a7e14dcfSSatish Balay {
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;
63a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
649566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n"));
65a7e14dcfSSatish Balay   }
66a7e14dcfSSatish Balay 
67a7e14dcfSSatish Balay   /* Initialized variables */
68a7e14dcfSSatish Balay   pert = nlsP->sval;
69a7e14dcfSSatish Balay 
702d9aa51bSJason Sarich   /* Number of times ksp stopped because of these reasons */
71a7e14dcfSSatish Balay   nlsP->ksp_atol = 0;
72a7e14dcfSSatish Balay   nlsP->ksp_rtol = 0;
73a7e14dcfSSatish Balay   nlsP->ksp_dtol = 0;
74a7e14dcfSSatish Balay   nlsP->ksp_ctol = 0;
75a7e14dcfSSatish Balay   nlsP->ksp_negc = 0;
76a7e14dcfSSatish Balay   nlsP->ksp_iter = 0;
77a7e14dcfSSatish Balay   nlsP->ksp_othr = 0;
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay   /* Initialize trust-region radius when using nash, stcg, or gltr
80ba7fe8fbSTodd Munson      Command automatically ignored for other methods
81ba7fe8fbSTodd Munson      Will be reset during the first iteration
82ba7fe8fbSTodd Munson   */
839566063dSJacob Faibussowitsch   PetscCall(KSPGetType(tao->ksp,&ksp_type));
849566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type,KSPNASH,&is_nash));
859566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type,KSPSTCG,&is_stcg));
869566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(ksp_type,KSPGLTR,&is_gltr));
871daac38eSTodd Munson 
889566063dSJacob Faibussowitsch   PetscCall(KSPCGSetRadius(tao->ksp,nlsP->max_radius));
89a7e14dcfSSatish Balay 
901daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
913c859ba3SBarry Smith     PetscCheck(tao->trust0 >= 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_OUTOFRANGE,"Initial radius negative");
92a7e14dcfSSatish Balay     tao->trust = tao->trust0;
93a7e14dcfSSatish Balay     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
94a7e14dcfSSatish Balay     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
95a7e14dcfSSatish Balay   }
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay   /* Check convergence criteria */
989566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
999566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
1003c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
101a7e14dcfSSatish Balay 
1023ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1039566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
1049566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
105*dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
1063ecd9318SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
107a7e14dcfSSatish Balay 
1080c51296cSAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
1099566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
1120c51296cSAlp Dener   if (is_bfgs) {
1130c51296cSAlp Dener     nlsP->bfgs_pre = pc;
1149566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(nlsP->bfgs_pre, &nlsP->M));
1159566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
1169566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
1179566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(nlsP->M, n, n, N, N));
1189566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(nlsP->M, tao->solution, tao->gradient));
1199566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(nlsP->M, &sym_set, &is_symmetric));
1203c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
1211baa6e33SBarry Smith   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc,PETSC_TRUE));
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay   /* Initialize trust-region radius.  The initialization is only performed
124a7e14dcfSSatish Balay      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
1251daac38eSTodd Munson   if (is_nash || is_stcg || is_gltr) {
126a7e14dcfSSatish Balay     switch (nlsP->init_type) {
127a7e14dcfSSatish Balay     case NLS_INIT_CONSTANT:
128a7e14dcfSSatish Balay       /* Use the initial radius specified */
129a7e14dcfSSatish Balay       break;
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay     case NLS_INIT_INTERPOLATION:
132a7e14dcfSSatish Balay       /* Use the initial radius specified */
133a7e14dcfSSatish Balay       max_radius = 0.0;
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay       for (j = 0; j < j_max; ++j) {
136a7e14dcfSSatish Balay         fmin = f;
137a7e14dcfSSatish Balay         sigma = 0.0;
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay         if (needH) {
1409566063dSJacob Faibussowitsch           PetscCall(TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre));
141a7e14dcfSSatish Balay           needH = 0;
142a7e14dcfSSatish Balay         }
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay         for (i = 0; i < i_max; ++i) {
1459566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution,nlsP->W));
1469566063dSJacob Faibussowitsch           PetscCall(VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient));
1479566063dSJacob Faibussowitsch           PetscCall(TaoComputeObjective(tao, nlsP->W, &ftrial));
148a7e14dcfSSatish Balay           if (PetscIsInfOrNanReal(ftrial)) {
149a7e14dcfSSatish Balay             tau = nlsP->gamma1_i;
15087f595a5SBarry Smith           } else {
151a7e14dcfSSatish Balay             if (ftrial < fmin) {
152a7e14dcfSSatish Balay               fmin = ftrial;
153a7e14dcfSSatish Balay               sigma = -tao->trust / gnorm;
154a7e14dcfSSatish Balay             }
155a7e14dcfSSatish Balay 
1569566063dSJacob Faibussowitsch             PetscCall(MatMult(tao->hessian, tao->gradient, nlsP->D));
1579566063dSJacob Faibussowitsch             PetscCall(VecDot(tao->gradient, nlsP->D, &prered));
158a7e14dcfSSatish Balay 
159a7e14dcfSSatish Balay             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
160a7e14dcfSSatish Balay             actred = f - ftrial;
16187f595a5SBarry Smith             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
162a7e14dcfSSatish Balay               kappa = 1.0;
16387f595a5SBarry Smith             } else {
164a7e14dcfSSatish Balay               kappa = actred / prered;
165a7e14dcfSSatish Balay             }
166a7e14dcfSSatish Balay 
167a7e14dcfSSatish Balay             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
168a7e14dcfSSatish Balay             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
169a7e14dcfSSatish Balay             tau_min = PetscMin(tau_1, tau_2);
170a7e14dcfSSatish Balay             tau_max = PetscMax(tau_1, tau_2);
171a7e14dcfSSatish Balay 
17218cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu1_i) {
173a7e14dcfSSatish Balay               /* Great agreement */
174a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
175a7e14dcfSSatish Balay 
176a7e14dcfSSatish Balay               if (tau_max < 1.0) {
177a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
17887f595a5SBarry Smith               } else if (tau_max > nlsP->gamma4_i) {
179a7e14dcfSSatish Balay                 tau = nlsP->gamma4_i;
18087f595a5SBarry Smith               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
181a7e14dcfSSatish Balay                 tau = tau_1;
18287f595a5SBarry Smith               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
183a7e14dcfSSatish Balay                 tau = tau_2;
18487f595a5SBarry Smith               } else {
185a7e14dcfSSatish Balay                 tau = tau_max;
186a7e14dcfSSatish Balay               }
18718cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= nlsP->mu2_i) {
188a7e14dcfSSatish Balay               /* Good agreement */
189a7e14dcfSSatish Balay               max_radius = PetscMax(max_radius, tao->trust);
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay               if (tau_max < nlsP->gamma2_i) {
192a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
19387f595a5SBarry Smith               } else if (tau_max > nlsP->gamma3_i) {
194a7e14dcfSSatish Balay                 tau = nlsP->gamma3_i;
19587f595a5SBarry Smith               } else {
196a7e14dcfSSatish Balay                 tau = tau_max;
197a7e14dcfSSatish Balay               }
19887f595a5SBarry Smith             } else {
199a7e14dcfSSatish Balay               /* Not good agreement */
200a7e14dcfSSatish Balay               if (tau_min > 1.0) {
201a7e14dcfSSatish Balay                 tau = nlsP->gamma2_i;
20287f595a5SBarry Smith               } else if (tau_max < nlsP->gamma1_i) {
203a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
20487f595a5SBarry Smith               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
205a7e14dcfSSatish Balay                 tau = nlsP->gamma1_i;
20687f595a5SBarry Smith               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
207a7e14dcfSSatish Balay                 tau = tau_1;
20887f595a5SBarry Smith               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
209a7e14dcfSSatish Balay                 tau = tau_2;
21087f595a5SBarry Smith               } else {
211a7e14dcfSSatish Balay                 tau = tau_max;
212a7e14dcfSSatish Balay               }
213a7e14dcfSSatish Balay             }
214a7e14dcfSSatish Balay           }
215a7e14dcfSSatish Balay           tao->trust = tau * tao->trust;
216a7e14dcfSSatish Balay         }
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay         if (fmin < f) {
219a7e14dcfSSatish Balay           f = fmin;
2209566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution,sigma,tao->gradient));
2219566063dSJacob Faibussowitsch           PetscCall(TaoComputeGradient(tao,tao->solution,tao->gradient));
222a7e14dcfSSatish Balay 
2239566063dSJacob Faibussowitsch           PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
2243c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute gradient generated Inf or NaN");
225a7e14dcfSSatish Balay           needH = 1;
226a7e14dcfSSatish Balay 
2279566063dSJacob Faibussowitsch           PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
2289566063dSJacob Faibussowitsch           PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
229*dbbe0bcdSBarry Smith           PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
2303ecd9318SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
231a7e14dcfSSatish Balay         }
232a7e14dcfSSatish Balay       }
233a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, max_radius);
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay       /* Modify the radius if it is too large or small */
236a7e14dcfSSatish Balay       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
237a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
238a7e14dcfSSatish Balay       break;
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay     default:
241a7e14dcfSSatish Balay       /* Norm of the first direction will initialize radius */
242a7e14dcfSSatish Balay       tao->trust = 0.0;
243a7e14dcfSSatish Balay       break;
244a7e14dcfSSatish Balay     }
245a7e14dcfSSatish Balay   }
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps*/
248a7e14dcfSSatish Balay   nlsP->newt = 0;
249a7e14dcfSSatish Balay   nlsP->bfgs = 0;
250a7e14dcfSSatish Balay   nlsP->grad = 0;
251a7e14dcfSSatish Balay 
252a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
2533ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
254e1e80dc8SAlp Dener     /* Call general purpose update function */
255*dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao,update, tao->niter, tao->user_update);
2568931d482SJason Sarich     ++tao->niter;
257ae93cb3cSJason Sarich     tao->ksp_its = 0;
258a7e14dcfSSatish Balay 
259a7e14dcfSSatish Balay     /* Compute the Hessian */
2601baa6e33SBarry Smith     if (needH) PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre));
261a7e14dcfSSatish Balay 
262a7e14dcfSSatish Balay     /* Shift the Hessian matrix */
263a7e14dcfSSatish Balay     if (pert > 0) {
2649566063dSJacob Faibussowitsch       PetscCall(MatShift(tao->hessian, pert));
265a7e14dcfSSatish Balay       if (tao->hessian != tao->hessian_pre) {
2669566063dSJacob Faibussowitsch         PetscCall(MatShift(tao->hessian_pre, pert));
267a7e14dcfSSatish Balay       }
268a7e14dcfSSatish Balay     }
269a7e14dcfSSatish Balay 
2700c51296cSAlp Dener     if (nlsP->bfgs_pre) {
2719566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
272a7e14dcfSSatish Balay       ++bfgsUpdates;
273a7e14dcfSSatish Balay     }
274a7e14dcfSSatish Balay 
275a7e14dcfSSatish Balay     /* Solve the Newton system of equations */
2769566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre));
2771daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
2789566063dSJacob Faibussowitsch       PetscCall(KSPCGSetRadius(tao->ksp,nlsP->max_radius));
2799566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
2809566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp,&kspits));
281a7e14dcfSSatish Balay       tao->ksp_its += kspits;
282ae93cb3cSJason Sarich       tao->ksp_tot_its += kspits;
2839566063dSJacob Faibussowitsch       PetscCall(KSPCGGetNormD(tao->ksp,&norm_d));
284a7e14dcfSSatish Balay 
285a7e14dcfSSatish Balay       if (0.0 == tao->trust) {
286a7e14dcfSSatish Balay         /* Radius was uninitialized; use the norm of the direction */
287a7e14dcfSSatish Balay         if (norm_d > 0.0) {
288a7e14dcfSSatish Balay           tao->trust = norm_d;
289a7e14dcfSSatish Balay 
290a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
291a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
292a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
29387f595a5SBarry Smith         } else {
294a7e14dcfSSatish Balay           /* The direction was bad; set radius to default value and re-solve
295a7e14dcfSSatish Balay              the trust-region subproblem to get a direction */
296a7e14dcfSSatish Balay           tao->trust = tao->trust0;
297a7e14dcfSSatish Balay 
298a7e14dcfSSatish Balay           /* Modify the radius if it is too large or small */
299a7e14dcfSSatish Balay           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
300a7e14dcfSSatish Balay           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
301a7e14dcfSSatish Balay 
3029566063dSJacob Faibussowitsch           PetscCall(KSPCGSetRadius(tao->ksp,nlsP->max_radius));
3039566063dSJacob Faibussowitsch           PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
3049566063dSJacob Faibussowitsch           PetscCall(KSPGetIterationNumber(tao->ksp,&kspits));
305a7e14dcfSSatish Balay           tao->ksp_its += kspits;
306ae93cb3cSJason Sarich           tao->ksp_tot_its += kspits;
3079566063dSJacob Faibussowitsch           PetscCall(KSPCGGetNormD(tao->ksp,&norm_d));
308a7e14dcfSSatish Balay 
3093c859ba3SBarry Smith           PetscCheck(norm_d != 0.0,PETSC_COMM_SELF,PETSC_ERR_PLIB, "Initial direction zero");
310a7e14dcfSSatish Balay         }
311a7e14dcfSSatish Balay       }
31287f595a5SBarry Smith     } else {
3139566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tao->gradient, nlsP->D));
3149566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
315a7e14dcfSSatish Balay       tao->ksp_its += kspits;
316ae93cb3cSJason Sarich       tao->ksp_tot_its += kspits;
317a7e14dcfSSatish Balay     }
3189566063dSJacob Faibussowitsch     PetscCall(VecScale(nlsP->D, -1.0));
3199566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(tao->ksp, &ksp_reason));
3200c51296cSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (nlsP->bfgs_pre)) {
321a7e14dcfSSatish Balay       /* Preconditioner is numerically indefinite; reset the
322a7e14dcfSSatish Balay          approximate if using BFGS preconditioning. */
3239566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
3249566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
325a7e14dcfSSatish Balay       bfgsUpdates = 1;
326a7e14dcfSSatish Balay     }
327a7e14dcfSSatish Balay 
328a7e14dcfSSatish Balay     if (KSP_CONVERGED_ATOL == ksp_reason) {
329a7e14dcfSSatish Balay       ++nlsP->ksp_atol;
33087f595a5SBarry Smith     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
331a7e14dcfSSatish Balay       ++nlsP->ksp_rtol;
33287f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
333a7e14dcfSSatish Balay       ++nlsP->ksp_ctol;
33487f595a5SBarry Smith     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
335a7e14dcfSSatish Balay       ++nlsP->ksp_negc;
33687f595a5SBarry Smith     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
337a7e14dcfSSatish Balay       ++nlsP->ksp_dtol;
33887f595a5SBarry Smith     } else if (KSP_DIVERGED_ITS == ksp_reason) {
339a7e14dcfSSatish Balay       ++nlsP->ksp_iter;
34087f595a5SBarry Smith     } else {
341a7e14dcfSSatish Balay       ++nlsP->ksp_othr;
342a7e14dcfSSatish Balay     }
343a7e14dcfSSatish Balay 
344a7e14dcfSSatish Balay     /* Check for success (descent direction) */
3459566063dSJacob Faibussowitsch     PetscCall(VecDot(nlsP->D, tao->gradient, &gdx));
346a7e14dcfSSatish Balay     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
347a7e14dcfSSatish Balay       /* Newton step is not descent or direction produced Inf or NaN
348a7e14dcfSSatish Balay          Update the perturbation for next time */
349a7e14dcfSSatish Balay       if (pert <= 0.0) {
350a7e14dcfSSatish Balay         /* Initialize the perturbation */
351a7e14dcfSSatish Balay         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
3521daac38eSTodd Munson         if (is_gltr) {
3539566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min));
354a7e14dcfSSatish Balay           pert = PetscMax(pert, -e_min);
355a7e14dcfSSatish Balay         }
35687f595a5SBarry Smith       } else {
357a7e14dcfSSatish Balay         /* Increase the perturbation */
358a7e14dcfSSatish Balay         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
359a7e14dcfSSatish Balay       }
360a7e14dcfSSatish Balay 
3610c51296cSAlp Dener       if (!nlsP->bfgs_pre) {
362a7e14dcfSSatish Balay         /* We don't have the bfgs matrix around and updated
363a7e14dcfSSatish Balay            Must use gradient direction in this case */
3649566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->gradient, nlsP->D));
3659566063dSJacob Faibussowitsch         PetscCall(VecScale(nlsP->D, -1.0));
366a7e14dcfSSatish Balay         ++nlsP->grad;
367a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
36887f595a5SBarry Smith       } else {
369a7e14dcfSSatish Balay         /* Attempt to use the BFGS direction */
3709566063dSJacob Faibussowitsch         PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
3719566063dSJacob Faibussowitsch         PetscCall(VecScale(nlsP->D, -1.0));
372a7e14dcfSSatish Balay 
373a7e14dcfSSatish Balay         /* Check for success (descent direction) */
3749566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, nlsP->D, &gdx));
375a7e14dcfSSatish Balay         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
376a7e14dcfSSatish Balay           /* BFGS direction is not descent or direction produced not a number
377a7e14dcfSSatish Balay              We can assert bfgsUpdates > 1 in this case because
378a7e14dcfSSatish Balay              the first solve produces the scaled gradient direction,
379a7e14dcfSSatish Balay              which is guaranteed to be descent */
380a7e14dcfSSatish Balay 
381a7e14dcfSSatish Balay           /* Use steepest descent direction (scaled) */
3829566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
3839566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
3849566063dSJacob Faibussowitsch           PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
3859566063dSJacob Faibussowitsch           PetscCall(VecScale(nlsP->D, -1.0));
386a7e14dcfSSatish Balay 
387a7e14dcfSSatish Balay           bfgsUpdates = 1;
3880c51296cSAlp Dener           ++nlsP->grad;
3890c51296cSAlp Dener           stepType = NLS_GRADIENT;
39087f595a5SBarry Smith         } else {
3919566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(nlsP->M, &bfgsUpdates));
392a7e14dcfSSatish Balay           if (1 == bfgsUpdates) {
393a7e14dcfSSatish Balay             /* The first BFGS direction is always the scaled gradient */
3940c51296cSAlp Dener             ++nlsP->grad;
3950c51296cSAlp Dener             stepType = NLS_GRADIENT;
39687f595a5SBarry Smith           } else {
397a7e14dcfSSatish Balay             ++nlsP->bfgs;
398a7e14dcfSSatish Balay             stepType = NLS_BFGS;
399a7e14dcfSSatish Balay           }
400a7e14dcfSSatish Balay         }
401a7e14dcfSSatish Balay       }
40287f595a5SBarry Smith     } else {
403a7e14dcfSSatish Balay       /* Computed Newton step is descent */
404a7e14dcfSSatish Balay       switch (ksp_reason) {
405a7e14dcfSSatish Balay       case KSP_DIVERGED_NANORINF:
406a7e14dcfSSatish Balay       case KSP_DIVERGED_BREAKDOWN:
407a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_MAT:
408a7e14dcfSSatish Balay       case KSP_DIVERGED_INDEFINITE_PC:
409a7e14dcfSSatish Balay       case KSP_CONVERGED_CG_NEG_CURVE:
410a7e14dcfSSatish Balay         /* Matrix or preconditioner is indefinite; increase perturbation */
411a7e14dcfSSatish Balay         if (pert <= 0.0) {
412a7e14dcfSSatish Balay           /* Initialize the perturbation */
413a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4141daac38eSTodd Munson           if (is_gltr) {
4159566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
416a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
417a7e14dcfSSatish Balay           }
41887f595a5SBarry Smith         } else {
419a7e14dcfSSatish Balay           /* Increase the perturbation */
420a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
421a7e14dcfSSatish Balay         }
422a7e14dcfSSatish Balay         break;
423a7e14dcfSSatish Balay 
424a7e14dcfSSatish Balay       default:
425a7e14dcfSSatish Balay         /* Newton step computation is good; decrease perturbation */
426a7e14dcfSSatish Balay         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
427a7e14dcfSSatish Balay         if (pert < nlsP->pmin) {
428a7e14dcfSSatish Balay           pert = 0.0;
429a7e14dcfSSatish Balay         }
430a7e14dcfSSatish Balay         break;
431a7e14dcfSSatish Balay       }
432a7e14dcfSSatish Balay 
433a7e14dcfSSatish Balay       ++nlsP->newt;
434a7e14dcfSSatish Balay       stepType = NLS_NEWTON;
435a7e14dcfSSatish Balay     }
436a7e14dcfSSatish Balay 
437a7e14dcfSSatish Balay     /* Perform the linesearch */
438a7e14dcfSSatish Balay     fold = f;
4399566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, nlsP->Xold));
4409566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, nlsP->Gold));
441a7e14dcfSSatish Balay 
4429566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason));
4439566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
444a7e14dcfSSatish Balay 
44587f595a5SBarry Smith     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
446a7e14dcfSSatish Balay       f = fold;
4479566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Xold, tao->solution));
4489566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Gold, tao->gradient));
449a7e14dcfSSatish Balay 
450a7e14dcfSSatish Balay       switch(stepType) {
451a7e14dcfSSatish Balay       case NLS_NEWTON:
452a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with Newton 1step
453a7e14dcfSSatish Balay            Update the perturbation for next time */
454a7e14dcfSSatish Balay         if (pert <= 0.0) {
455a7e14dcfSSatish Balay           /* Initialize the perturbation */
456a7e14dcfSSatish Balay           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
4571daac38eSTodd Munson           if (is_gltr) {
4589566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min));
459a7e14dcfSSatish Balay             pert = PetscMax(pert, -e_min);
460a7e14dcfSSatish Balay           }
46187f595a5SBarry Smith         } else {
462a7e14dcfSSatish Balay           /* Increase the perturbation */
463a7e14dcfSSatish Balay           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
464a7e14dcfSSatish Balay         }
465a7e14dcfSSatish Balay 
4660c51296cSAlp Dener         if (!nlsP->bfgs_pre) {
467a7e14dcfSSatish Balay           /* We don't have the bfgs matrix around and being updated
468a7e14dcfSSatish Balay              Must use gradient direction in this case */
4699566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->gradient, nlsP->D));
470a7e14dcfSSatish Balay           ++nlsP->grad;
471a7e14dcfSSatish Balay           stepType = NLS_GRADIENT;
47287f595a5SBarry Smith         } else {
473a7e14dcfSSatish Balay           /* Attempt to use the BFGS direction */
4749566063dSJacob Faibussowitsch           PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
475a7e14dcfSSatish Balay           /* Check for success (descent direction) */
4769566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->solution, nlsP->D, &gdx));
477a7e14dcfSSatish Balay           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
478a7e14dcfSSatish Balay             /* BFGS direction is not descent or direction produced not a number
479a7e14dcfSSatish Balay                We can assert bfgsUpdates > 1 in this case
480a7e14dcfSSatish Balay                Use steepest descent direction (scaled) */
4819566063dSJacob Faibussowitsch             PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
4829566063dSJacob Faibussowitsch             PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
4839566063dSJacob Faibussowitsch             PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
484a7e14dcfSSatish Balay 
485a7e14dcfSSatish Balay             bfgsUpdates = 1;
4860c51296cSAlp Dener             ++nlsP->grad;
4870c51296cSAlp Dener             stepType = NLS_GRADIENT;
48887f595a5SBarry Smith           } else {
489a7e14dcfSSatish Balay             if (1 == bfgsUpdates) {
490a7e14dcfSSatish Balay               /* The first BFGS direction is always the scaled gradient */
4910c51296cSAlp Dener               ++nlsP->grad;
4920c51296cSAlp Dener               stepType = NLS_GRADIENT;
49387f595a5SBarry Smith             } else {
494a7e14dcfSSatish Balay               ++nlsP->bfgs;
495a7e14dcfSSatish Balay               stepType = NLS_BFGS;
496a7e14dcfSSatish Balay             }
497a7e14dcfSSatish Balay           }
498a7e14dcfSSatish Balay         }
499a7e14dcfSSatish Balay         break;
500a7e14dcfSSatish Balay 
501a7e14dcfSSatish Balay       case NLS_BFGS:
502a7e14dcfSSatish Balay         /* Can only enter if pc_type == NLS_PC_BFGS
503a7e14dcfSSatish Balay            Failed to obtain acceptable iterate with BFGS step
504a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
5059566063dSJacob Faibussowitsch         PetscCall(MatLMVMReset(nlsP->M, PETSC_FALSE));
5069566063dSJacob Faibussowitsch         PetscCall(MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient));
5079566063dSJacob Faibussowitsch         PetscCall(MatSolve(nlsP->M, tao->gradient, nlsP->D));
508a7e14dcfSSatish Balay 
509a7e14dcfSSatish Balay         bfgsUpdates = 1;
510a7e14dcfSSatish Balay         ++nlsP->grad;
511a7e14dcfSSatish Balay         stepType = NLS_GRADIENT;
512a7e14dcfSSatish Balay         break;
513a7e14dcfSSatish Balay       }
5149566063dSJacob Faibussowitsch       PetscCall(VecScale(nlsP->D, -1.0));
515a7e14dcfSSatish Balay 
5169566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason));
5179566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full));
5189566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
519a7e14dcfSSatish Balay     }
520a7e14dcfSSatish Balay 
52187f595a5SBarry Smith     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
522a7e14dcfSSatish Balay       /* Failed to find an improving point */
523a7e14dcfSSatish Balay       f = fold;
5249566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Xold, tao->solution));
5259566063dSJacob Faibussowitsch       PetscCall(VecCopy(nlsP->Gold, tao->gradient));
526a7e14dcfSSatish Balay       step = 0.0;
527a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
528a7e14dcfSSatish Balay       break;
529a7e14dcfSSatish Balay     }
530a7e14dcfSSatish Balay 
531a7e14dcfSSatish Balay     /* Update trust region radius */
5321daac38eSTodd Munson     if (is_nash || is_stcg || is_gltr) {
533a7e14dcfSSatish Balay       switch(nlsP->update_type) {
534a7e14dcfSSatish Balay       case NLS_UPDATE_STEP:
535a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
536a7e14dcfSSatish Balay           if (step < nlsP->nu1) {
537a7e14dcfSSatish Balay             /* Very bad step taken; reduce radius */
538a7e14dcfSSatish Balay             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
53987f595a5SBarry Smith           } else if (step < nlsP->nu2) {
540a7e14dcfSSatish Balay             /* Reasonably bad step taken; reduce radius */
541a7e14dcfSSatish Balay             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
54287f595a5SBarry Smith           } else if (step < nlsP->nu3) {
543a7e14dcfSSatish Balay             /*  Reasonable step was taken; leave radius alone */
544a7e14dcfSSatish Balay             if (nlsP->omega3 < 1.0) {
545a7e14dcfSSatish Balay               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
54687f595a5SBarry Smith             } else if (nlsP->omega3 > 1.0) {
547a7e14dcfSSatish Balay               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
548a7e14dcfSSatish Balay             }
54987f595a5SBarry Smith           } else if (step < nlsP->nu4) {
550a7e14dcfSSatish Balay             /*  Full step taken; increase the radius */
551a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
55287f595a5SBarry Smith           } else {
553a7e14dcfSSatish Balay             /*  More than full step taken; increase the radius */
554a7e14dcfSSatish Balay             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
555a7e14dcfSSatish Balay           }
55687f595a5SBarry Smith         } else {
557a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
558a7e14dcfSSatish Balay           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
559a7e14dcfSSatish Balay         }
560a7e14dcfSSatish Balay         break;
561a7e14dcfSSatish Balay 
562a7e14dcfSSatish Balay       case NLS_UPDATE_REDUCTION:
563a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
564a7e14dcfSSatish Balay           /*  Get predicted reduction */
565a7e14dcfSSatish Balay 
5669566063dSJacob Faibussowitsch           PetscCall(KSPCGGetObjFcn(tao->ksp,&prered));
567a7e14dcfSSatish Balay           if (prered >= 0.0) {
568a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
569a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
570a7e14dcfSSatish Balay             /*  be rejected! */
571a7e14dcfSSatish Balay             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
57287f595a5SBarry Smith           } else {
573a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
574a7e14dcfSSatish Balay               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
57587f595a5SBarry Smith             } else {
576a7e14dcfSSatish Balay               /*  Compute and actual reduction */
577a7e14dcfSSatish Balay               actred = fold - f_full;
578a7e14dcfSSatish Balay               prered = -prered;
579a7e14dcfSSatish Balay               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
580a7e14dcfSSatish Balay                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
581a7e14dcfSSatish Balay                 kappa = 1.0;
58287f595a5SBarry Smith               } else {
583a7e14dcfSSatish Balay                 kappa = actred / prered;
584a7e14dcfSSatish Balay               }
585a7e14dcfSSatish Balay 
586a7e14dcfSSatish Balay               /*  Accept of reject the step and update radius */
587a7e14dcfSSatish Balay               if (kappa < nlsP->eta1) {
588a7e14dcfSSatish Balay                 /*  Very bad step */
589a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
59087f595a5SBarry Smith               } else if (kappa < nlsP->eta2) {
591a7e14dcfSSatish Balay                 /*  Marginal bad step */
592a7e14dcfSSatish Balay                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
59387f595a5SBarry Smith               } else if (kappa < nlsP->eta3) {
594a7e14dcfSSatish Balay                 /*  Reasonable step */
595a7e14dcfSSatish Balay                 if (nlsP->alpha3 < 1.0) {
596a7e14dcfSSatish Balay                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
59787f595a5SBarry Smith                 } else if (nlsP->alpha3 > 1.0) {
598a7e14dcfSSatish Balay                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
599a7e14dcfSSatish Balay                 }
60087f595a5SBarry Smith               } else if (kappa < nlsP->eta4) {
601a7e14dcfSSatish Balay                 /*  Good step */
602a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
60387f595a5SBarry Smith               } else {
604a7e14dcfSSatish Balay                 /*  Very good step */
605a7e14dcfSSatish Balay                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
606a7e14dcfSSatish Balay               }
607a7e14dcfSSatish Balay             }
608a7e14dcfSSatish Balay           }
60987f595a5SBarry Smith         } else {
610a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
611a7e14dcfSSatish Balay           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
612a7e14dcfSSatish Balay         }
613a7e14dcfSSatish Balay         break;
614a7e14dcfSSatish Balay 
615a7e14dcfSSatish Balay       default:
616a7e14dcfSSatish Balay         if (stepType == NLS_NEWTON) {
6179566063dSJacob Faibussowitsch           PetscCall(KSPCGGetObjFcn(tao->ksp,&prered));
618a7e14dcfSSatish Balay           if (prered >= 0.0) {
619a7e14dcfSSatish Balay             /*  The predicted reduction has the wrong sign.  This cannot */
620a7e14dcfSSatish Balay             /*  happen in infinite precision arithmetic.  Step should */
621a7e14dcfSSatish Balay             /*  be rejected! */
622a7e14dcfSSatish Balay             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
62387f595a5SBarry Smith           } else {
624a7e14dcfSSatish Balay             if (PetscIsInfOrNanReal(f_full)) {
625a7e14dcfSSatish Balay               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
62687f595a5SBarry Smith             } else {
627a7e14dcfSSatish Balay               actred = fold - f_full;
628a7e14dcfSSatish Balay               prered = -prered;
62987f595a5SBarry Smith               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
630a7e14dcfSSatish Balay                 kappa = 1.0;
63187f595a5SBarry Smith               } else {
632a7e14dcfSSatish Balay                 kappa = actred / prered;
633a7e14dcfSSatish Balay               }
634a7e14dcfSSatish Balay 
635a7e14dcfSSatish Balay               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
636a7e14dcfSSatish Balay               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
637a7e14dcfSSatish Balay               tau_min = PetscMin(tau_1, tau_2);
638a7e14dcfSSatish Balay               tau_max = PetscMax(tau_1, tau_2);
639a7e14dcfSSatish Balay 
640a7e14dcfSSatish Balay               if (kappa >= 1.0 - nlsP->mu1) {
641a7e14dcfSSatish Balay                 /*  Great agreement */
642a7e14dcfSSatish Balay                 if (tau_max < 1.0) {
643a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
64487f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma4) {
645a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
64687f595a5SBarry Smith                 } else {
647a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
648a7e14dcfSSatish Balay                 }
64987f595a5SBarry Smith               } else if (kappa >= 1.0 - nlsP->mu2) {
650a7e14dcfSSatish Balay                 /*  Good agreement */
651a7e14dcfSSatish Balay 
652a7e14dcfSSatish Balay                 if (tau_max < nlsP->gamma2) {
653a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
65487f595a5SBarry Smith                 } else if (tau_max > nlsP->gamma3) {
655a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
65687f595a5SBarry Smith                 } else if (tau_max < 1.0) {
657a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
65887f595a5SBarry Smith                 } else {
659a7e14dcfSSatish Balay                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
660a7e14dcfSSatish Balay                 }
66187f595a5SBarry Smith               } else {
662a7e14dcfSSatish Balay                 /*  Not good agreement */
663a7e14dcfSSatish Balay                 if (tau_min > 1.0) {
664a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
66587f595a5SBarry Smith                 } else if (tau_max < nlsP->gamma1) {
666a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
66787f595a5SBarry Smith                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
668a7e14dcfSSatish Balay                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
66987f595a5SBarry Smith                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
670a7e14dcfSSatish Balay                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
67187f595a5SBarry Smith                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
672a7e14dcfSSatish Balay                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
67387f595a5SBarry Smith                 } else {
674a7e14dcfSSatish Balay                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
675a7e14dcfSSatish Balay                 }
676a7e14dcfSSatish Balay               }
677a7e14dcfSSatish Balay             }
678a7e14dcfSSatish Balay           }
67987f595a5SBarry Smith         } else {
680a7e14dcfSSatish Balay           /*  Newton step was not good; reduce the radius */
681a7e14dcfSSatish Balay           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
682a7e14dcfSSatish Balay         }
683a7e14dcfSSatish Balay         break;
684a7e14dcfSSatish Balay       }
685a7e14dcfSSatish Balay 
686a7e14dcfSSatish Balay       /*  The radius may have been increased; modify if it is too large */
687a7e14dcfSSatish Balay       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
688a7e14dcfSSatish Balay     }
689a7e14dcfSSatish Balay 
690a7e14dcfSSatish Balay     /*  Check for termination */
6919566063dSJacob Faibussowitsch     PetscCall(TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm));
6923c859ba3SBarry Smith     PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER,"User provided compute function generated Not-a-Number");
693a7e14dcfSSatish Balay     needH = 1;
6949566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its));
6959566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao,tao->niter,f,gnorm,0.0,step));
696*dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao,convergencetest ,tao->cnvP);
697a7e14dcfSSatish Balay   }
698a7e14dcfSSatish Balay   PetscFunctionReturn(0);
699a7e14dcfSSatish Balay }
700a7e14dcfSSatish Balay 
701a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
702441846f8SBarry Smith static PetscErrorCode TaoSetUp_NLS(Tao tao)
703a7e14dcfSSatish Balay {
704a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
705a7e14dcfSSatish Balay 
706a7e14dcfSSatish Balay   PetscFunctionBegin;
7079566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution,&tao->gradient));
7089566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
7099566063dSJacob Faibussowitsch   if (!nlsP->W) PetscCall(VecDuplicate(tao->solution,&nlsP->W));
7109566063dSJacob Faibussowitsch   if (!nlsP->D) PetscCall(VecDuplicate(tao->solution,&nlsP->D));
7119566063dSJacob Faibussowitsch   if (!nlsP->Xold) PetscCall(VecDuplicate(tao->solution,&nlsP->Xold));
7129566063dSJacob Faibussowitsch   if (!nlsP->Gold) PetscCall(VecDuplicate(tao->solution,&nlsP->Gold));
71383c8fe1dSLisandro Dalcin   nlsP->M = NULL;
71483c8fe1dSLisandro Dalcin   nlsP->bfgs_pre = NULL;
715a7e14dcfSSatish Balay   PetscFunctionReturn(0);
716a7e14dcfSSatish Balay }
717a7e14dcfSSatish Balay 
718a7e14dcfSSatish Balay /*------------------------------------------------------------*/
719441846f8SBarry Smith static PetscErrorCode TaoDestroy_NLS(Tao tao)
720a7e14dcfSSatish Balay {
721a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
722a7e14dcfSSatish Balay 
723a7e14dcfSSatish Balay   PetscFunctionBegin;
724a7e14dcfSSatish Balay   if (tao->setupcalled) {
7259566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->D));
7269566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->W));
7279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->Xold));
7289566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nlsP->Gold));
729a7e14dcfSSatish Balay   }
730a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
7319566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
732a7e14dcfSSatish Balay   PetscFunctionReturn(0);
733a7e14dcfSSatish Balay }
734a7e14dcfSSatish Balay 
735a7e14dcfSSatish Balay /*------------------------------------------------------------*/
736*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_NLS(Tao tao,PetscOptionItems *PetscOptionsObject)
737a7e14dcfSSatish Balay {
738a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
739a7e14dcfSSatish Balay 
740a7e14dcfSSatish Balay   PetscFunctionBegin;
741d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Newton line search method for unconstrained optimization");
7429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, NULL));
7439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, NULL));
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL));
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL));
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL));
7479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL));
7489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL));
7499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL));
7509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL));
7519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL));
7529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL));
7539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL));
7549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL));
7559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL));
7569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL));
7579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL));
7589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL));
7599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL));
7609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL));
7619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL));
7629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL));
7639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL));
7649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL));
7659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL));
7669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL));
7679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL));
7689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL));
7699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL));
7709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL));
7719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL));
7729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL));
7739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL));
7749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL));
7759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL));
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL));
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL));
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL));
7799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL));
7809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL));
7819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL));
7829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL));
7839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL));
7849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL));
7859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL));
7869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL));
7879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL));
7889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL));
789d0609cedSBarry Smith   PetscOptionsHeadEnd();
7909566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
7919566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
792a7e14dcfSSatish Balay   PetscFunctionReturn(0);
793a7e14dcfSSatish Balay }
794a7e14dcfSSatish Balay 
795a7e14dcfSSatish Balay /*------------------------------------------------------------*/
796441846f8SBarry Smith static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
797a7e14dcfSSatish Balay {
798a7e14dcfSSatish Balay   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
799a7e14dcfSSatish Balay   PetscBool      isascii;
800a7e14dcfSSatish Balay 
801a7e14dcfSSatish Balay   PetscFunctionBegin;
8029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
803a7e14dcfSSatish Balay   if (isascii) {
8049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
80563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", nlsP->newt));
80663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", nlsP->bfgs));
80763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", nlsP->grad));
808a7e14dcfSSatish Balay 
80963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp atol: %" PetscInt_FMT "\n", nlsP->ksp_atol));
81063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %" PetscInt_FMT "\n", nlsP->ksp_rtol));
81163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %" PetscInt_FMT "\n", nlsP->ksp_ctol));
81263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp negc: %" PetscInt_FMT "\n", nlsP->ksp_negc));
81363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %" PetscInt_FMT "\n", nlsP->ksp_dtol));
81463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp iter: %" PetscInt_FMT "\n", nlsP->ksp_iter));
81563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "nls ksp othr: %" PetscInt_FMT "\n", nlsP->ksp_othr));
8169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
817a7e14dcfSSatish Balay   }
818a7e14dcfSSatish Balay   PetscFunctionReturn(0);
819a7e14dcfSSatish Balay }
820a7e14dcfSSatish Balay 
821a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
8224aa34175SJason Sarich /*MC
8234aa34175SJason Sarich   TAONLS - Newton's method with linesearch for unconstrained minimization.
8244aa34175SJason Sarich   At each iteration, the Newton line search method solves the symmetric
8254aa34175SJason Sarich   system of equations to obtain the step diretion dk:
8264aa34175SJason Sarich               Hk dk = -gk
8274aa34175SJason Sarich   a More-Thuente line search is applied on the direction dk to approximately
8284aa34175SJason Sarich   solve
8294aa34175SJason Sarich               min_t f(xk + t d_k)
8304aa34175SJason Sarich 
8314aa34175SJason Sarich     Options Database Keys:
8329d0a60b2SAlp Dener + -tao_nls_init_type - "constant","direction","interpolation"
8334aa34175SJason Sarich . -tao_nls_update_type - "step","direction","interpolation"
8344aa34175SJason Sarich . -tao_nls_sval - perturbation starting value
8354aa34175SJason Sarich . -tao_nls_imin - minimum initial perturbation
8364aa34175SJason Sarich . -tao_nls_imax - maximum initial perturbation
8374aa34175SJason Sarich . -tao_nls_pmin - minimum perturbation
8384aa34175SJason Sarich . -tao_nls_pmax - maximum perturbation
8394aa34175SJason Sarich . -tao_nls_pgfac - growth factor
8404aa34175SJason Sarich . -tao_nls_psfac - shrink factor
8414aa34175SJason Sarich . -tao_nls_imfac - initial merit factor
8424aa34175SJason Sarich . -tao_nls_pmgfac - merit growth factor
8434aa34175SJason Sarich . -tao_nls_pmsfac - merit shrink factor
8444aa34175SJason Sarich . -tao_nls_eta1 - poor steplength; reduce radius
8454aa34175SJason Sarich . -tao_nls_eta2 - reasonable steplength; leave radius
8464aa34175SJason Sarich . -tao_nls_eta3 - good steplength; increase readius
8474aa34175SJason Sarich . -tao_nls_eta4 - excellent steplength; greatly increase radius
8484aa34175SJason Sarich . -tao_nls_alpha1 - alpha1 reduction
8494aa34175SJason Sarich . -tao_nls_alpha2 - alpha2 reduction
8504aa34175SJason Sarich . -tao_nls_alpha3 - alpha3 reduction
8514aa34175SJason Sarich . -tao_nls_alpha4 - alpha4 reduction
8524aa34175SJason Sarich . -tao_nls_alpha - alpha5 reduction
8534aa34175SJason Sarich . -tao_nls_mu1 - mu1 interpolation update
8544aa34175SJason Sarich . -tao_nls_mu2 - mu2 interpolation update
8554aa34175SJason Sarich . -tao_nls_gamma1 - gamma1 interpolation update
8564aa34175SJason Sarich . -tao_nls_gamma2 - gamma2 interpolation update
8574aa34175SJason Sarich . -tao_nls_gamma3 - gamma3 interpolation update
8584aa34175SJason Sarich . -tao_nls_gamma4 - gamma4 interpolation update
8594aa34175SJason Sarich . -tao_nls_theta - theta interpolation update
8604aa34175SJason Sarich . -tao_nls_omega1 - omega1 step update
8614aa34175SJason Sarich . -tao_nls_omega2 - omega2 step update
8624aa34175SJason Sarich . -tao_nls_omega3 - omega3 step update
8634aa34175SJason Sarich . -tao_nls_omega4 - omega4 step update
8644aa34175SJason Sarich . -tao_nls_omega5 - omega5 step update
8651522df2eSJason Sarich . -tao_nls_mu1_i -  mu1 interpolation init factor
8661522df2eSJason Sarich . -tao_nls_mu2_i -  mu2 interpolation init factor
8671522df2eSJason Sarich . -tao_nls_gamma1_i -  gamma1 interpolation init factor
8681522df2eSJason Sarich . -tao_nls_gamma2_i -  gamma2 interpolation init factor
8691522df2eSJason Sarich . -tao_nls_gamma3_i -  gamma3 interpolation init factor
8701522df2eSJason Sarich . -tao_nls_gamma4_i -  gamma4 interpolation init factor
8711522df2eSJason Sarich - -tao_nls_theta_i -  theta interpolation init factor
8721eb8069cSJason Sarich 
8731eb8069cSJason Sarich   Level: beginner
8741522df2eSJason Sarich M*/
8754aa34175SJason Sarich 
876728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao)
877a7e14dcfSSatish Balay {
878a7e14dcfSSatish Balay   TAO_NLS        *nlsP;
8798caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
880a7e14dcfSSatish Balay 
881a7e14dcfSSatish Balay   PetscFunctionBegin;
8829566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&nlsP));
883a7e14dcfSSatish Balay 
884a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_NLS;
885a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_NLS;
886a7e14dcfSSatish Balay   tao->ops->view = TaoView_NLS;
887a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
888a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_NLS;
889a7e14dcfSSatish Balay 
8906552cf8aSJason Sarich   /* Override default settings (unless already changed) */
8916552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
8926552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 100.0;
8936552cf8aSJason Sarich 
894a7e14dcfSSatish Balay   tao->data = (void*)nlsP;
895a7e14dcfSSatish Balay 
896a7e14dcfSSatish Balay   nlsP->sval   = 0.0;
897a7e14dcfSSatish Balay   nlsP->imin   = 1.0e-4;
898a7e14dcfSSatish Balay   nlsP->imax   = 1.0e+2;
899a7e14dcfSSatish Balay   nlsP->imfac  = 1.0e-1;
900a7e14dcfSSatish Balay 
901a7e14dcfSSatish Balay   nlsP->pmin   = 1.0e-12;
902a7e14dcfSSatish Balay   nlsP->pmax   = 1.0e+2;
903a7e14dcfSSatish Balay   nlsP->pgfac  = 1.0e+1;
904a7e14dcfSSatish Balay   nlsP->psfac  = 4.0e-1;
905a7e14dcfSSatish Balay   nlsP->pmgfac = 1.0e-1;
906a7e14dcfSSatish Balay   nlsP->pmsfac = 1.0e-1;
907a7e14dcfSSatish Balay 
908a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on steplength */
909a7e14dcfSSatish Balay   nlsP->nu1 = 0.25;
910a7e14dcfSSatish Balay   nlsP->nu2 = 0.50;
911a7e14dcfSSatish Balay   nlsP->nu3 = 1.00;
912a7e14dcfSSatish Balay   nlsP->nu4 = 1.25;
913a7e14dcfSSatish Balay 
914a7e14dcfSSatish Balay   nlsP->omega1 = 0.25;
915a7e14dcfSSatish Balay   nlsP->omega2 = 0.50;
916a7e14dcfSSatish Balay   nlsP->omega3 = 1.00;
917a7e14dcfSSatish Balay   nlsP->omega4 = 2.00;
918a7e14dcfSSatish Balay   nlsP->omega5 = 4.00;
919a7e14dcfSSatish Balay 
920a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on reduction */
921a7e14dcfSSatish Balay   nlsP->eta1 = 1.0e-4;
922a7e14dcfSSatish Balay   nlsP->eta2 = 0.25;
923a7e14dcfSSatish Balay   nlsP->eta3 = 0.50;
924a7e14dcfSSatish Balay   nlsP->eta4 = 0.90;
925a7e14dcfSSatish Balay 
926a7e14dcfSSatish Balay   nlsP->alpha1 = 0.25;
927a7e14dcfSSatish Balay   nlsP->alpha2 = 0.50;
928a7e14dcfSSatish Balay   nlsP->alpha3 = 1.00;
929a7e14dcfSSatish Balay   nlsP->alpha4 = 2.00;
930a7e14dcfSSatish Balay   nlsP->alpha5 = 4.00;
931a7e14dcfSSatish Balay 
932a7e14dcfSSatish Balay   /*  Default values for trust-region radius update based on interpolation */
933a7e14dcfSSatish Balay   nlsP->mu1 = 0.10;
934a7e14dcfSSatish Balay   nlsP->mu2 = 0.50;
935a7e14dcfSSatish Balay 
936a7e14dcfSSatish Balay   nlsP->gamma1 = 0.25;
937a7e14dcfSSatish Balay   nlsP->gamma2 = 0.50;
938a7e14dcfSSatish Balay   nlsP->gamma3 = 2.00;
939a7e14dcfSSatish Balay   nlsP->gamma4 = 4.00;
940a7e14dcfSSatish Balay 
941a7e14dcfSSatish Balay   nlsP->theta = 0.05;
942a7e14dcfSSatish Balay 
943a7e14dcfSSatish Balay   /*  Default values for trust region initialization based on interpolation */
944a7e14dcfSSatish Balay   nlsP->mu1_i = 0.35;
945a7e14dcfSSatish Balay   nlsP->mu2_i = 0.50;
946a7e14dcfSSatish Balay 
947a7e14dcfSSatish Balay   nlsP->gamma1_i = 0.0625;
948a7e14dcfSSatish Balay   nlsP->gamma2_i = 0.5;
949a7e14dcfSSatish Balay   nlsP->gamma3_i = 2.0;
950a7e14dcfSSatish Balay   nlsP->gamma4_i = 5.0;
951a7e14dcfSSatish Balay 
952a7e14dcfSSatish Balay   nlsP->theta_i = 0.25;
953a7e14dcfSSatish Balay 
954a7e14dcfSSatish Balay   /*  Remaining parameters */
955a7e14dcfSSatish Balay   nlsP->min_radius = 1.0e-10;
956a7e14dcfSSatish Balay   nlsP->max_radius = 1.0e10;
957a7e14dcfSSatish Balay   nlsP->epsilon = 1.0e-6;
958a7e14dcfSSatish Balay 
959a7e14dcfSSatish Balay   nlsP->init_type       = NLS_INIT_INTERPOLATION;
960a7e14dcfSSatish Balay   nlsP->update_type     = NLS_UPDATE_STEP;
961a7e14dcfSSatish Balay 
9629566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch));
9639566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch,(PetscObject)tao,1));
9649566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type));
9659566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao));
9669566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix));
967a7e14dcfSSatish Balay 
968a7e14dcfSSatish Balay   /*  Set linear solver to default for symmetric matrices */
9699566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp));
9709566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1));
9719566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix));
9729566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp,"tao_nls_"));
9739566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp,KSPSTCG));
974a7e14dcfSSatish Balay   PetscFunctionReturn(0);
975a7e14dcfSSatish Balay }
976