xref: /petsc/src/tao/unconstrained/impls/nls/nls.c (revision 6aeb920f53dc9577f1df4ff00bd447ebf24efec0)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/matrix/lmvmmat.h>
3 #include <../src/tao/unconstrained/impls/nls/nls.h>
4 
5 #include <petscksp.h>
6 #include <petscpc.h>
7 #include <petsc-private/kspimpl.h>
8 #include <petsc-private/pcimpl.h>
9 
10 #define NLS_KSP_CG      0
11 #define NLS_KSP_NASH    1
12 #define NLS_KSP_STCG    2
13 #define NLS_KSP_GLTR    3
14 #define NLS_KSP_PETSC   4
15 #define NLS_KSP_TYPES   5
16 
17 #define NLS_PC_NONE     0
18 #define NLS_PC_AHESS    1
19 #define NLS_PC_BFGS     2
20 #define NLS_PC_PETSC    3
21 #define NLS_PC_TYPES    4
22 
23 #define BFGS_SCALE_AHESS        0
24 #define BFGS_SCALE_PHESS        1
25 #define BFGS_SCALE_BFGS         2
26 #define BFGS_SCALE_TYPES        3
27 
28 #define NLS_INIT_CONSTANT         0
29 #define NLS_INIT_DIRECTION        1
30 #define NLS_INIT_INTERPOLATION    2
31 #define NLS_INIT_TYPES            3
32 
33 #define NLS_UPDATE_STEP           0
34 #define NLS_UPDATE_REDUCTION      1
35 #define NLS_UPDATE_INTERPOLATION  2
36 #define NLS_UPDATE_TYPES          3
37 
38 static const char *NLS_KSP[64] = {"cg", "nash", "stcg", "gltr", "petsc"};
39 
40 static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"};
41 
42 static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
43 
44 static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"};
45 
46 static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"};
47 
48 static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x);
49 /* Routine for BFGS preconditioner
50 
51 
52  Implements Newton's Method with a line search approach for solving
53  unconstrained minimization problems.  A More'-Thuente line search
54  is used to guarantee that the bfgs preconditioner remains positive
55  definite.
56 
57  The method can shift the Hessian matrix.  The shifting procedure is
58  adapted from the PATH algorithm for solving complementarity
59  problems.
60 
61  The linear system solve should be done with a conjugate gradient
62  method, although any method can be used. */
63 
64 #define NLS_NEWTON              0
65 #define NLS_BFGS                1
66 #define NLS_SCALED_GRADIENT     2
67 #define NLS_GRADIENT            3
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "TaoSolve_NLS"
71 static PetscErrorCode TaoSolve_NLS(Tao tao)
72 {
73   PetscErrorCode               ierr;
74   TAO_NLS                      *nlsP = (TAO_NLS *)tao->data;
75   PC                           pc;
76   KSPConvergedReason           ksp_reason;
77   TaoLineSearchConvergedReason ls_reason;
78   TaoConvergedReason           reason;
79 
80   PetscReal                    fmin, ftrial, f_full, prered, actred, kappa, sigma;
81   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
82   PetscReal                    f, fold, gdx, gnorm, pert;
83   PetscReal                    step = 1.0;
84   PetscReal                    delta;
85   PetscReal                    norm_d = 0.0, e_min;
86 
87   PetscInt                     stepType;
88   PetscInt                     iter = 0;
89   PetscInt                     bfgsUpdates = 0;
90   PetscInt                     n,N,kspits;
91   PetscInt                     needH;
92 
93   PetscInt                     i_max = 5;
94   PetscInt                     j_max = 1;
95   PetscInt                     i, j;
96 
97   PetscFunctionBegin;
98   if (tao->XL || tao->XU || tao->ops->computebounds) {
99     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr);
100   }
101 
102   /* Initialized variables */
103   pert = nlsP->sval;
104 
105   nlsP->ksp_atol = 0;
106   nlsP->ksp_rtol = 0;
107   nlsP->ksp_dtol = 0;
108   nlsP->ksp_ctol = 0;
109   nlsP->ksp_negc = 0;
110   nlsP->ksp_iter = 0;
111   nlsP->ksp_othr = 0;
112 
113   /* Modify the linear solver to a trust region method if desired */
114   switch(nlsP->ksp_type) {
115   case NLS_KSP_CG:
116     ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr);
117     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
118     break;
119 
120   case NLS_KSP_NASH:
121     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
122     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
123     break;
124 
125   case NLS_KSP_STCG:
126     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
127     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
128     break;
129 
130   case NLS_KSP_GLTR:
131     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
132     ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
133     break;
134 
135   default:
136     /* Use the method set by the ksp_type */
137     break;
138   }
139 
140   /* Initialize trust-region radius when using nash, stcg, or gltr
141    Will be reset during the first iteration */
142   if (NLS_KSP_NASH == nlsP->ksp_type) {
143     ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
144   } else if (NLS_KSP_STCG == nlsP->ksp_type) {
145     ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
146   } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
147     ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
148   }
149 
150   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
151     tao->trust = tao->trust0;
152 
153     if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative");
154 
155     /* Modify the radius if it is too large or small */
156     tao->trust = PetscMax(tao->trust, nlsP->min_radius);
157     tao->trust = PetscMin(tao->trust, nlsP->max_radius);
158   }
159 
160   /* Get vectors we will need */
161   if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) {
162     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
163     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
164     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr);
165     ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr);
166   }
167 
168   /* Check convergence criteria */
169   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
170   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
171   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
172   needH = 1;
173 
174   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
175   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
176 
177   /* create vectors for the limited memory preconditioner */
178   if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) {
179     if (!nlsP->Diag) {
180       ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr);
181     }
182   }
183 
184   /* Modify the preconditioner to use the bfgs approximation */
185   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
186   switch(nlsP->pc_type) {
187   case NLS_PC_NONE:
188     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
189     if (pc->ops->setfromoptions) {
190       (*pc->ops->setfromoptions)(pc);
191     }
192     break;
193 
194   case NLS_PC_AHESS:
195     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
196     if (pc->ops->setfromoptions) {
197       (*pc->ops->setfromoptions)(pc);
198     }
199     ierr = PCJacobiSetUseAbs(pc);CHKERRQ(ierr);
200     break;
201 
202   case NLS_PC_BFGS:
203     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
204     if (pc->ops->setfromoptions) {
205       (*pc->ops->setfromoptions)(pc);
206     }
207     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
208     ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr);
209     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
210     break;
211 
212   default:
213     /* Use the pc method set by pc_type */
214     break;
215   }
216 
217   /* Initialize trust-region radius.  The initialization is only performed
218      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
219   if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
220     switch(nlsP->init_type) {
221     case NLS_INIT_CONSTANT:
222       /* Use the initial radius specified */
223       break;
224 
225     case NLS_INIT_INTERPOLATION:
226       /* Use the initial radius specified */
227       max_radius = 0.0;
228 
229       for (j = 0; j < j_max; ++j) {
230         fmin = f;
231         sigma = 0.0;
232 
233         if (needH) {
234           ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
235           needH = 0;
236         }
237 
238         for (i = 0; i < i_max; ++i) {
239           ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr);
240           ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr);
241           ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr);
242           if (PetscIsInfOrNanReal(ftrial)) {
243             tau = nlsP->gamma1_i;
244           } else {
245             if (ftrial < fmin) {
246               fmin = ftrial;
247               sigma = -tao->trust / gnorm;
248             }
249 
250             ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr);
251             ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr);
252 
253             prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
254             actred = f - ftrial;
255             if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
256               kappa = 1.0;
257             } else {
258               kappa = actred / prered;
259             }
260 
261             tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred);
262             tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred);
263             tau_min = PetscMin(tau_1, tau_2);
264             tau_max = PetscMax(tau_1, tau_2);
265 
266             if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) {
267               /* Great agreement */
268               max_radius = PetscMax(max_radius, tao->trust);
269 
270               if (tau_max < 1.0) {
271                 tau = nlsP->gamma3_i;
272               } else if (tau_max > nlsP->gamma4_i) {
273                 tau = nlsP->gamma4_i;
274               } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) {
275                 tau = tau_1;
276               } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) {
277                 tau = tau_2;
278               } else {
279                 tau = tau_max;
280               }
281             } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) {
282               /* Good agreement */
283               max_radius = PetscMax(max_radius, tao->trust);
284 
285               if (tau_max < nlsP->gamma2_i) {
286                 tau = nlsP->gamma2_i;
287               } else if (tau_max > nlsP->gamma3_i) {
288                 tau = nlsP->gamma3_i;
289               } else {
290                 tau = tau_max;
291               }
292             } else {
293               /* Not good agreement */
294               if (tau_min > 1.0) {
295                 tau = nlsP->gamma2_i;
296               } else if (tau_max < nlsP->gamma1_i) {
297                 tau = nlsP->gamma1_i;
298               } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) {
299                 tau = nlsP->gamma1_i;
300               } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
301                 tau = tau_1;
302               } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) {
303                 tau = tau_2;
304               } else {
305                 tau = tau_max;
306               }
307             }
308           }
309           tao->trust = tau * tao->trust;
310         }
311 
312         if (fmin < f) {
313           f = fmin;
314           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
315           ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr);
316 
317           ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
318           if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
319           needH = 1;
320 
321           ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
322           if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
323         }
324       }
325       tao->trust = PetscMax(tao->trust, max_radius);
326 
327       /* Modify the radius if it is too large or small */
328       tao->trust = PetscMax(tao->trust, nlsP->min_radius);
329       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
330       break;
331 
332     default:
333       /* Norm of the first direction will initialize radius */
334       tao->trust = 0.0;
335       break;
336     }
337   }
338 
339   /* Set initial scaling for the BFGS preconditioner
340      This step is done after computing the initial trust-region radius
341      since the function value may have decreased */
342   if (NLS_PC_BFGS == nlsP->pc_type) {
343     if (f != 0.0) {
344       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
345     } else {
346       delta = 2.0 / (gnorm*gnorm);
347     }
348     ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
349   }
350 
351   /* Set counter for gradient/reset steps*/
352   nlsP->newt = 0;
353   nlsP->bfgs = 0;
354   nlsP->sgrad = 0;
355   nlsP->grad = 0;
356 
357   /* Have not converged; continue with Newton method */
358   while (reason == TAO_CONTINUE_ITERATING) {
359     ++iter;
360     tao->ksp_its=0;
361 
362     /* Compute the Hessian */
363     if (needH) {
364       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
365       needH = 0;
366     }
367 
368     if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) {
369       /* Obtain diagonal for the bfgs preconditioner  */
370       ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
371       ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
372       ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
373       ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
374     }
375 
376     /* Shift the Hessian matrix */
377     if (pert > 0) {
378       ierr = MatShift(tao->hessian, pert);
379       if (tao->hessian != tao->hessian_pre) {
380         ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr);
381       }
382     }
383 
384     if (NLS_PC_BFGS == nlsP->pc_type) {
385       if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) {
386         /* Obtain diagonal for the bfgs preconditioner  */
387         ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr);
388         ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr);
389         ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr);
390         ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr);
391       }
392       /* Update the limited memory preconditioner */
393       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
394       ++bfgsUpdates;
395     }
396 
397     /* Solve the Newton system of equations */
398     ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
399     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type ||  NLS_KSP_GLTR == nlsP->ksp_type) {
400 
401       if (NLS_KSP_NASH == nlsP->ksp_type) {
402         ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
403       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
404          ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
405       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
406         ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
407       }
408 
409       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
410       ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
411       tao->ksp_its+=kspits;
412       tao->ksp_tot_its+=kspits;
413 
414       if (NLS_KSP_NASH == nlsP->ksp_type) {
415         ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
416       } else if (NLS_KSP_STCG == nlsP->ksp_type) {
417          ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
418       } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
419         ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
420       }
421 
422       if (0.0 == tao->trust) {
423         /* Radius was uninitialized; use the norm of the direction */
424         if (norm_d > 0.0) {
425           tao->trust = norm_d;
426 
427           /* Modify the radius if it is too large or small */
428           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
429           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
430         } else {
431           /* The direction was bad; set radius to default value and re-solve
432              the trust-region subproblem to get a direction */
433           tao->trust = tao->trust0;
434 
435           /* Modify the radius if it is too large or small */
436           tao->trust = PetscMax(tao->trust, nlsP->min_radius);
437           tao->trust = PetscMin(tao->trust, nlsP->max_radius);
438 
439           if (NLS_KSP_NASH == nlsP->ksp_type) {
440             ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
441           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
442             ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
443           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
444             ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr);
445           }
446 
447           ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
448           ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
449           tao->ksp_its+=kspits;
450           tao->ksp_tot_its+=kspits;
451           if (NLS_KSP_NASH == nlsP->ksp_type) {
452             ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
453           } else if (NLS_KSP_STCG == nlsP->ksp_type) {
454             ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
455           } else if (NLS_KSP_GLTR == nlsP->ksp_type) {
456             ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr);
457           }
458 
459           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
460         }
461       }
462     } else {
463       ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr);
464       ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
465       tao->ksp_its += kspits;
466       tao->ksp_tot_its+=kspits;
467     }
468     ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
469     ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
470     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) {
471       /* Preconditioner is numerically indefinite; reset the
472          approximate if using BFGS preconditioning. */
473 
474       if (f != 0.0) {
475         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
476       } else {
477         delta = 2.0 / (gnorm*gnorm);
478       }
479       ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr);
480       ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
481       ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
482       bfgsUpdates = 1;
483     }
484 
485     if (KSP_CONVERGED_ATOL == ksp_reason) {
486       ++nlsP->ksp_atol;
487     } else if (KSP_CONVERGED_RTOL == ksp_reason) {
488       ++nlsP->ksp_rtol;
489     } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
490       ++nlsP->ksp_ctol;
491     } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
492       ++nlsP->ksp_negc;
493     } else if (KSP_DIVERGED_DTOL == ksp_reason) {
494       ++nlsP->ksp_dtol;
495     } else if (KSP_DIVERGED_ITS == ksp_reason) {
496       ++nlsP->ksp_iter;
497     } else {
498       ++nlsP->ksp_othr;
499     }
500 
501     /* Check for success (descent direction) */
502     ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr);
503     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
504       /* Newton step is not descent or direction produced Inf or NaN
505          Update the perturbation for next time */
506       if (pert <= 0.0) {
507         /* Initialize the perturbation */
508         pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
509         if (NLS_KSP_GLTR == nlsP->ksp_type) {
510           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
511           pert = PetscMax(pert, -e_min);
512         }
513       } else {
514         /* Increase the perturbation */
515         pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
516       }
517 
518       if (NLS_PC_BFGS != nlsP->pc_type) {
519         /* We don't have the bfgs matrix around and updated
520            Must use gradient direction in this case */
521         ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
522         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
523         ++nlsP->grad;
524         stepType = NLS_GRADIENT;
525       } else {
526         /* Attempt to use the BFGS direction */
527         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
528         ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
529 
530         /* Check for success (descent direction) */
531         ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr);
532         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
533           /* BFGS direction is not descent or direction produced not a number
534              We can assert bfgsUpdates > 1 in this case because
535              the first solve produces the scaled gradient direction,
536              which is guaranteed to be descent */
537 
538           /* Use steepest descent direction (scaled) */
539 
540           if (f != 0.0) {
541             delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
542           } else {
543             delta = 2.0 / (gnorm*gnorm);
544           }
545           ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
546           ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
547           ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
548           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
549           ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
550 
551           bfgsUpdates = 1;
552           ++nlsP->sgrad;
553           stepType = NLS_SCALED_GRADIENT;
554         } else {
555           if (1 == bfgsUpdates) {
556             /* The first BFGS direction is always the scaled gradient */
557             ++nlsP->sgrad;
558             stepType = NLS_SCALED_GRADIENT;
559           } else {
560             ++nlsP->bfgs;
561             stepType = NLS_BFGS;
562           }
563         }
564       }
565     } else {
566       /* Computed Newton step is descent */
567       switch (ksp_reason) {
568       case KSP_DIVERGED_NANORINF:
569       case KSP_DIVERGED_BREAKDOWN:
570       case KSP_DIVERGED_INDEFINITE_MAT:
571       case KSP_DIVERGED_INDEFINITE_PC:
572       case KSP_CONVERGED_CG_NEG_CURVE:
573         /* Matrix or preconditioner is indefinite; increase perturbation */
574         if (pert <= 0.0) {
575           /* Initialize the perturbation */
576           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
577           if (NLS_KSP_GLTR == nlsP->ksp_type) {
578             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
579             pert = PetscMax(pert, -e_min);
580           }
581         } else {
582           /* Increase the perturbation */
583           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
584         }
585         break;
586 
587       default:
588         /* Newton step computation is good; decrease perturbation */
589         pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm);
590         if (pert < nlsP->pmin) {
591           pert = 0.0;
592         }
593         break;
594       }
595 
596       ++nlsP->newt;
597       stepType = NLS_NEWTON;
598     }
599 
600     /* Perform the linesearch */
601     fold = f;
602     ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr);
603     ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr);
604 
605     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
606     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
607 
608     while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) {      /* Linesearch failed */
609       f = fold;
610       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
611       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
612 
613       switch(stepType) {
614       case NLS_NEWTON:
615         /* Failed to obtain acceptable iterate with Newton 1step
616            Update the perturbation for next time */
617         if (pert <= 0.0) {
618           /* Initialize the perturbation */
619           pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm));
620           if (NLS_KSP_GLTR == nlsP->ksp_type) {
621             ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
622             pert = PetscMax(pert, -e_min);
623           }
624         } else {
625           /* Increase the perturbation */
626           pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm));
627         }
628 
629         if (NLS_PC_BFGS != nlsP->pc_type) {
630           /* We don't have the bfgs matrix around and being updated
631              Must use gradient direction in this case */
632           ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr);
633           ++nlsP->grad;
634           stepType = NLS_GRADIENT;
635         } else {
636           /* Attempt to use the BFGS direction */
637           ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
638           /* Check for success (descent direction) */
639           ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr);
640           if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
641             /* BFGS direction is not descent or direction produced not a number
642                We can assert bfgsUpdates > 1 in this case
643                Use steepest descent direction (scaled) */
644 
645             if (f != 0.0) {
646               delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
647             } else {
648               delta = 2.0 / (gnorm*gnorm);
649             }
650             ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
651             ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
652             ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
653             ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
654 
655             bfgsUpdates = 1;
656             ++nlsP->sgrad;
657             stepType = NLS_SCALED_GRADIENT;
658           } else {
659             if (1 == bfgsUpdates) {
660               /* The first BFGS direction is always the scaled gradient */
661               ++nlsP->sgrad;
662               stepType = NLS_SCALED_GRADIENT;
663             } else {
664               ++nlsP->bfgs;
665               stepType = NLS_BFGS;
666             }
667           }
668         }
669         break;
670 
671       case NLS_BFGS:
672         /* Can only enter if pc_type == NLS_PC_BFGS
673            Failed to obtain acceptable iterate with BFGS step
674            Attempt to use the scaled gradient direction */
675 
676         if (f != 0.0) {
677           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
678         } else {
679           delta = 2.0 / (gnorm*gnorm);
680         }
681         ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr);
682         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
683         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
684         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
685 
686         bfgsUpdates = 1;
687         ++nlsP->sgrad;
688         stepType = NLS_SCALED_GRADIENT;
689         break;
690 
691       case NLS_SCALED_GRADIENT:
692         /* Can only enter if pc_type == NLS_PC_BFGS
693            The scaled gradient step did not produce a new iterate;
694            attemp to use the gradient direction.
695            Need to make sure we are not using a different diagonal scaling */
696 
697         ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr);
698         ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr);
699         ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr);
700         ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
701         ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr);
702 
703         bfgsUpdates = 1;
704         ++nlsP->grad;
705         stepType = NLS_GRADIENT;
706         break;
707       }
708       ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr);
709 
710       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr);
711       ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr);
712       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
713     }
714 
715     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
716       /* Failed to find an improving point */
717       f = fold;
718       ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr);
719       ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr);
720       step = 0.0;
721       reason = TAO_DIVERGED_LS_FAILURE;
722       tao->reason = TAO_DIVERGED_LS_FAILURE;
723       break;
724     }
725 
726     /* Update trust region radius */
727     if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) {
728       switch(nlsP->update_type) {
729       case NLS_UPDATE_STEP:
730         if (stepType == NLS_NEWTON) {
731           if (step < nlsP->nu1) {
732             /* Very bad step taken; reduce radius */
733             tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
734           } else if (step < nlsP->nu2) {
735             /* Reasonably bad step taken; reduce radius */
736             tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust);
737           } else if (step < nlsP->nu3) {
738             /*  Reasonable step was taken; leave radius alone */
739             if (nlsP->omega3 < 1.0) {
740               tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust);
741             } else if (nlsP->omega3 > 1.0) {
742               tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust);
743             }
744           } else if (step < nlsP->nu4) {
745             /*  Full step taken; increase the radius */
746             tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust);
747           } else {
748             /*  More than full step taken; increase the radius */
749             tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust);
750           }
751         } else {
752           /*  Newton step was not good; reduce the radius */
753           tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust);
754         }
755         break;
756 
757       case NLS_UPDATE_REDUCTION:
758         if (stepType == NLS_NEWTON) {
759           /*  Get predicted reduction */
760 
761           if (NLS_KSP_STCG == nlsP->ksp_type) {
762             ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);
763           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
764             ierr = KSPNASHGetObjFcn(tao->ksp,&prered);
765           } else {
766             ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
767           }
768 
769           if (prered >= 0.0) {
770             /*  The predicted reduction has the wrong sign.  This cannot */
771             /*  happen in infinite precision arithmetic.  Step should */
772             /*  be rejected! */
773             tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
774           } else {
775             if (PetscIsInfOrNanReal(f_full)) {
776               tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
777             } else {
778               /*  Compute and actual reduction */
779               actred = fold - f_full;
780               prered = -prered;
781               if ((PetscAbsScalar(actred) <= nlsP->epsilon) &&
782                   (PetscAbsScalar(prered) <= nlsP->epsilon)) {
783                 kappa = 1.0;
784               } else {
785                 kappa = actred / prered;
786               }
787 
788               /*  Accept of reject the step and update radius */
789               if (kappa < nlsP->eta1) {
790                 /*  Very bad step */
791                 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d);
792               } else if (kappa < nlsP->eta2) {
793                 /*  Marginal bad step */
794                 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d);
795               } else if (kappa < nlsP->eta3) {
796                 /*  Reasonable step */
797                 if (nlsP->alpha3 < 1.0) {
798                   tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust);
799                 } else if (nlsP->alpha3 > 1.0) {
800                   tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust);
801                 }
802               } else if (kappa < nlsP->eta4) {
803                 /*  Good step */
804                 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust);
805               } else {
806                 /*  Very good step */
807                 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust);
808               }
809             }
810           }
811         } else {
812           /*  Newton step was not good; reduce the radius */
813           tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust);
814         }
815         break;
816 
817       default:
818         if (stepType == NLS_NEWTON) {
819 
820           if (NLS_KSP_STCG == nlsP->ksp_type) {
821               ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);
822           } else if (NLS_KSP_NASH == nlsP->ksp_type)  {
823               ierr = KSPNASHGetObjFcn(tao->ksp,&prered);
824           } else {
825               ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
826           }
827           if (prered >= 0.0) {
828             /*  The predicted reduction has the wrong sign.  This cannot */
829             /*  happen in infinite precision arithmetic.  Step should */
830             /*  be rejected! */
831             tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
832           } else {
833             if (PetscIsInfOrNanReal(f_full)) {
834               tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
835             } else {
836               actred = fold - f_full;
837               prered = -prered;
838               if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) {
839                 kappa = 1.0;
840               } else {
841                 kappa = actred / prered;
842               }
843 
844               tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred);
845               tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred);
846               tau_min = PetscMin(tau_1, tau_2);
847               tau_max = PetscMax(tau_1, tau_2);
848 
849               if (kappa >= 1.0 - nlsP->mu1) {
850                 /*  Great agreement */
851                 if (tau_max < 1.0) {
852                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
853                 } else if (tau_max > nlsP->gamma4) {
854                   tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d);
855                 } else {
856                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
857                 }
858               } else if (kappa >= 1.0 - nlsP->mu2) {
859                 /*  Good agreement */
860 
861                 if (tau_max < nlsP->gamma2) {
862                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
863                 } else if (tau_max > nlsP->gamma3) {
864                   tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d);
865                 } else if (tau_max < 1.0) {
866                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
867                 } else {
868                   tao->trust = PetscMax(tao->trust, tau_max * norm_d);
869                 }
870               } else {
871                 /*  Not good agreement */
872                 if (tau_min > 1.0) {
873                   tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d);
874                 } else if (tau_max < nlsP->gamma1) {
875                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
876                 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) {
877                   tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d);
878                 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) {
879                   tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
880                 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) {
881                   tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
882                 } else {
883                   tao->trust = tau_max * PetscMin(tao->trust, norm_d);
884                 }
885               }
886             }
887           }
888         } else {
889           /*  Newton step was not good; reduce the radius */
890           tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust);
891         }
892         break;
893       }
894 
895       /*  The radius may have been increased; modify if it is too large */
896       tao->trust = PetscMin(tao->trust, nlsP->max_radius);
897     }
898 
899     /*  Check for termination */
900     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
901     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
902     needH = 1;
903     ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
904   }
905   PetscFunctionReturn(0);
906 }
907 
908 /* ---------------------------------------------------------- */
909 #undef __FUNCT__
910 #define __FUNCT__ "TaoSetUp_NLS"
911 static PetscErrorCode TaoSetUp_NLS(Tao tao)
912 {
913   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
914   PetscErrorCode ierr;
915 
916   PetscFunctionBegin;
917   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
918   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
919   if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);}
920   if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);}
921   if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);}
922   if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);}
923   nlsP->Diag = 0;
924   nlsP->M = 0;
925   PetscFunctionReturn(0);
926 }
927 
928 /*------------------------------------------------------------*/
929 #undef __FUNCT__
930 #define __FUNCT__ "TaoDestroy_NLS"
931 static PetscErrorCode TaoDestroy_NLS(Tao tao)
932 {
933   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
934   PetscErrorCode ierr;
935 
936   PetscFunctionBegin;
937   if (tao->setupcalled) {
938     ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr);
939     ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr);
940     ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr);
941     ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr);
942   }
943   ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr);
944   ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr);
945   ierr = PetscFree(tao->data);CHKERRQ(ierr);
946   PetscFunctionReturn(0);
947 }
948 
949 /*------------------------------------------------------------*/
950 #undef __FUNCT__
951 #define __FUNCT__ "TaoSetFromOptions_NLS"
952 static PetscErrorCode TaoSetFromOptions_NLS(Tao tao)
953 {
954   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
955   PetscErrorCode ierr;
956 
957   PetscFunctionBegin;
958   ierr = PetscOptionsHead("Newton line search method for unconstrained optimization");CHKERRQ(ierr);
959   ierr = PetscOptionsEList("-tao_nls_ksp_type", "ksp type", "", NLS_KSP, NLS_KSP_TYPES, NLS_KSP[nlsP->ksp_type], &nlsP->ksp_type, 0);CHKERRQ(ierr);
960   ierr = PetscOptionsEList("-tao_nls_pc_type", "pc type", "", NLS_PC, NLS_PC_TYPES, NLS_PC[nlsP->pc_type], &nlsP->pc_type, 0);CHKERRQ(ierr);
961   ierr = PetscOptionsEList("-tao_nls_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[nlsP->bfgs_scale_type], &nlsP->bfgs_scale_type, 0);CHKERRQ(ierr);
962   ierr = PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, 0);CHKERRQ(ierr);
963   ierr = PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, 0);CHKERRQ(ierr);
964  ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval, 0);CHKERRQ(ierr);
965   ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin, 0);CHKERRQ(ierr);
966   ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax, 0);CHKERRQ(ierr);
967   ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac, 0);CHKERRQ(ierr);
968   ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin, 0);CHKERRQ(ierr);
969   ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax, 0);CHKERRQ(ierr);
970   ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac, 0);CHKERRQ(ierr);
971   ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac, 0);CHKERRQ(ierr);
972   ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac, 0);CHKERRQ(ierr);
973   ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac, 0);CHKERRQ(ierr);
974   ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1, 0);CHKERRQ(ierr);
975   ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2, 0);CHKERRQ(ierr);
976   ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3, 0);CHKERRQ(ierr);
977   ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4, 0);CHKERRQ(ierr);
978   ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1, 0);CHKERRQ(ierr);
979   ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2, 0);CHKERRQ(ierr);
980   ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3, 0);CHKERRQ(ierr);
981   ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4, 0);CHKERRQ(ierr);
982   ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5, 0);CHKERRQ(ierr);
983   ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1, 0);CHKERRQ(ierr);
984   ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2, 0);CHKERRQ(ierr);
985   ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3, 0);CHKERRQ(ierr);
986   ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4, 0);CHKERRQ(ierr);
987   ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1, 0);CHKERRQ(ierr);
988   ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2, 0);CHKERRQ(ierr);
989   ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3, 0);CHKERRQ(ierr);
990   ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4, 0);CHKERRQ(ierr);
991   ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5, 0);CHKERRQ(ierr);
992   ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i, 0);CHKERRQ(ierr);
993   ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i, 0);CHKERRQ(ierr);
994   ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i, 0);CHKERRQ(ierr);
995   ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i, 0);CHKERRQ(ierr);
996   ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i, 0);CHKERRQ(ierr);
997   ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i, 0);CHKERRQ(ierr);
998   ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i, 0);CHKERRQ(ierr);
999   ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1, 0);CHKERRQ(ierr);
1000   ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2, 0);CHKERRQ(ierr);
1001   ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1, 0);CHKERRQ(ierr);
1002   ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2, 0);CHKERRQ(ierr);
1003   ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3, 0);CHKERRQ(ierr);
1004   ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4, 0);CHKERRQ(ierr);
1005   ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta, 0);CHKERRQ(ierr);
1006   ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius, 0);CHKERRQ(ierr);
1007   ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius, 0);CHKERRQ(ierr);
1008   ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon, 0);CHKERRQ(ierr);
1009   ierr = PetscOptionsTail();CHKERRQ(ierr);
1010   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1011   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 
1016 /*------------------------------------------------------------*/
1017 #undef __FUNCT__
1018 #define __FUNCT__ "TaoView_NLS"
1019 static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer)
1020 {
1021   TAO_NLS        *nlsP = (TAO_NLS *)tao->data;
1022   PetscInt       nrejects;
1023   PetscBool      isascii;
1024   PetscErrorCode ierr;
1025 
1026   PetscFunctionBegin;
1027   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1028   if (isascii) {
1029     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1030     if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) {
1031       ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr);
1032       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1033     }
1034     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr);
1035     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr);
1036     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr);
1037     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr);
1038 
1039     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr);
1040     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr);
1041     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr);
1042     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr);
1043     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr);
1044     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr);
1045     ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr);
1046     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1047   }
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 /* ---------------------------------------------------------- */
1052 /*MC
1053   TAONLS - Newton's method with linesearch for unconstrained minimization.
1054   At each iteration, the Newton line search method solves the symmetric
1055   system of equations to obtain the step diretion dk:
1056               Hk dk = -gk
1057   a More-Thuente line search is applied on the direction dk to approximately
1058   solve
1059               min_t f(xk + t d_k)
1060 
1061     Options Database Keys:
1062 + -tao_nls_ksp_type - "cg","nash","stcg","gltr","petsc"
1063 . -tao_nls_pc_type - "none","ahess","bfgs","petsc"
1064 . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs"
1065 . -tao_nls_init_type - "constant","direction","interpolation"
1066 . -tao_nls_update_type - "step","direction","interpolation"
1067 . -tao_nls_sval - perturbation starting value
1068 . -tao_nls_imin - minimum initial perturbation
1069 . -tao_nls_imax - maximum initial perturbation
1070 . -tao_nls_pmin - minimum perturbation
1071 . -tao_nls_pmax - maximum perturbation
1072 . -tao_nls_pgfac - growth factor
1073 . -tao_nls_psfac - shrink factor
1074 . -tao_nls_imfac - initial merit factor
1075 . -tao_nls_pmgfac - merit growth factor
1076 . -tao_nls_pmsfac - merit shrink factor
1077 . -tao_nls_eta1 - poor steplength; reduce radius
1078 . -tao_nls_eta2 - reasonable steplength; leave radius
1079 . -tao_nls_eta3 - good steplength; increase readius
1080 . -tao_nls_eta4 - excellent steplength; greatly increase radius
1081 . -tao_nls_alpha1 - alpha1 reduction
1082 . -tao_nls_alpha2 - alpha2 reduction
1083 . -tao_nls_alpha3 - alpha3 reduction
1084 . -tao_nls_alpha4 - alpha4 reduction
1085 . -tao_nls_alpha - alpha5 reduction
1086 . -tao_nls_mu1 - mu1 interpolation update
1087 . -tao_nls_mu2 - mu2 interpolation update
1088 . -tao_nls_gamma1 - gamma1 interpolation update
1089 . -tao_nls_gamma2 - gamma2 interpolation update
1090 . -tao_nls_gamma3 - gamma3 interpolation update
1091 . -tao_nls_gamma4 - gamma4 interpolation update
1092 . -tao_nls_theta - theta interpolation update
1093 . -tao_nls_omega1 - omega1 step update
1094 . -tao_nls_omega2 - omega2 step update
1095 . -tao_nls_omega3 - omega3 step update
1096 . -tao_nls_omega4 - omega4 step update
1097 . -tao_nls_omega5 - omega5 step update
1098 . -tao_nls_mu1_i -  mu1 interpolation init factor
1099 . -tao_nls_mu2_i -  mu2 interpolation init factor
1100 . -tao_nls_gamma1_i -  gamma1 interpolation init factor
1101 . -tao_nls_gamma2_i -  gamma2 interpolation init factor
1102 . -tao_nls_gamma3_i -  gamma3 interpolation init factor
1103 . -tao_nls_gamma4_i -  gamma4 interpolation init factor
1104 - -tao_nls_theta_i -  theta interpolation init factor
1105 
1106   Level: beginner
1107 M*/
1108 
1109 EXTERN_C_BEGIN
1110 #undef __FUNCT__
1111 #define __FUNCT__ "TaoCreate_NLS"
1112 PetscErrorCode TaoCreate_NLS(Tao tao)
1113 {
1114   TAO_NLS        *nlsP;
1115   const char     *morethuente_type = TAOLINESEARCHMT;
1116   PetscErrorCode ierr;
1117 
1118   PetscFunctionBegin;
1119   ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr);
1120 
1121   tao->ops->setup = TaoSetUp_NLS;
1122   tao->ops->solve = TaoSolve_NLS;
1123   tao->ops->view = TaoView_NLS;
1124   tao->ops->setfromoptions = TaoSetFromOptions_NLS;
1125   tao->ops->destroy = TaoDestroy_NLS;
1126 
1127   tao->max_it = 50;
1128 #if defined(PETSC_USE_REAL_SINGLE)
1129   tao->fatol = 1e-5;
1130   tao->frtol = 1e-5;
1131 #else
1132   tao->fatol = 1e-10;
1133   tao->frtol = 1e-10;
1134 #endif
1135   tao->data = (void*)nlsP;
1136   tao->trust0 = 100.0;
1137 
1138   nlsP->sval   = 0.0;
1139   nlsP->imin   = 1.0e-4;
1140   nlsP->imax   = 1.0e+2;
1141   nlsP->imfac  = 1.0e-1;
1142 
1143   nlsP->pmin   = 1.0e-12;
1144   nlsP->pmax   = 1.0e+2;
1145   nlsP->pgfac  = 1.0e+1;
1146   nlsP->psfac  = 4.0e-1;
1147   nlsP->pmgfac = 1.0e-1;
1148   nlsP->pmsfac = 1.0e-1;
1149 
1150   /*  Default values for trust-region radius update based on steplength */
1151   nlsP->nu1 = 0.25;
1152   nlsP->nu2 = 0.50;
1153   nlsP->nu3 = 1.00;
1154   nlsP->nu4 = 1.25;
1155 
1156   nlsP->omega1 = 0.25;
1157   nlsP->omega2 = 0.50;
1158   nlsP->omega3 = 1.00;
1159   nlsP->omega4 = 2.00;
1160   nlsP->omega5 = 4.00;
1161 
1162   /*  Default values for trust-region radius update based on reduction */
1163   nlsP->eta1 = 1.0e-4;
1164   nlsP->eta2 = 0.25;
1165   nlsP->eta3 = 0.50;
1166   nlsP->eta4 = 0.90;
1167 
1168   nlsP->alpha1 = 0.25;
1169   nlsP->alpha2 = 0.50;
1170   nlsP->alpha3 = 1.00;
1171   nlsP->alpha4 = 2.00;
1172   nlsP->alpha5 = 4.00;
1173 
1174   /*  Default values for trust-region radius update based on interpolation */
1175   nlsP->mu1 = 0.10;
1176   nlsP->mu2 = 0.50;
1177 
1178   nlsP->gamma1 = 0.25;
1179   nlsP->gamma2 = 0.50;
1180   nlsP->gamma3 = 2.00;
1181   nlsP->gamma4 = 4.00;
1182 
1183   nlsP->theta = 0.05;
1184 
1185   /*  Default values for trust region initialization based on interpolation */
1186   nlsP->mu1_i = 0.35;
1187   nlsP->mu2_i = 0.50;
1188 
1189   nlsP->gamma1_i = 0.0625;
1190   nlsP->gamma2_i = 0.5;
1191   nlsP->gamma3_i = 2.0;
1192   nlsP->gamma4_i = 5.0;
1193 
1194   nlsP->theta_i = 0.25;
1195 
1196   /*  Remaining parameters */
1197   nlsP->min_radius = 1.0e-10;
1198   nlsP->max_radius = 1.0e10;
1199   nlsP->epsilon = 1.0e-6;
1200 
1201   nlsP->ksp_type        = NLS_KSP_STCG;
1202   nlsP->pc_type         = NLS_PC_BFGS;
1203   nlsP->bfgs_scale_type = BFGS_SCALE_PHESS;
1204   nlsP->init_type       = NLS_INIT_INTERPOLATION;
1205   nlsP->update_type     = NLS_UPDATE_STEP;
1206 
1207   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1208   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1209   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1210 
1211   /*  Set linear solver to default for symmetric matrices */
1212   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1213   PetscFunctionReturn(0);
1214 }
1215 EXTERN_C_END
1216 
1217 #undef __FUNCT__
1218 #define __FUNCT__ "MatLMVMSolveShell"
1219 static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
1220 {
1221   PetscErrorCode ierr;
1222   Mat            M;
1223 
1224   PetscFunctionBegin;
1225   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1226   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
1227   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
1228   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
1229   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
1230   PetscFunctionReturn(0);
1231 }
1232