xref: /petsc/src/tao/unconstrained/impls/ntr/ntr.c (revision 3e7ff0edd3573be01c8c0fa32db97c3db8fa5c8d)
1 #include <../src/tao/matrix/lmvmmat.h>
2 #include <../src/tao/unconstrained/impls/ntr/ntr.h>
3 
4 #include <petscksp.h>
5 #include <petscpc.h>
6 #include <petsc-private/kspimpl.h>
7 #include <petsc-private/pcimpl.h>
8 
9 #define NTR_KSP_NASH    0
10 #define NTR_KSP_STCG    1
11 #define NTR_KSP_GLTR    2
12 #define NTR_KSP_TYPES   3
13 
14 #define NTR_PC_NONE     0
15 #define NTR_PC_AHESS    1
16 #define NTR_PC_BFGS     2
17 #define NTR_PC_PETSC    3
18 #define NTR_PC_TYPES    4
19 
20 #define BFGS_SCALE_AHESS   0
21 #define BFGS_SCALE_BFGS    1
22 #define BFGS_SCALE_TYPES   2
23 
24 #define NTR_INIT_CONSTANT         0
25 #define NTR_INIT_DIRECTION        1
26 #define NTR_INIT_INTERPOLATION    2
27 #define NTR_INIT_TYPES            3
28 
29 #define NTR_UPDATE_REDUCTION      0
30 #define NTR_UPDATE_INTERPOLATION  1
31 #define NTR_UPDATE_TYPES          2
32 
33 static const char *NTR_KSP[64] = {  "nash", "stcg", "gltr"};
34 
35 static const char *NTR_PC[64] = {  "none", "ahess", "bfgs", "petsc"};
36 
37 static const char *BFGS_SCALE[64] = {  "ahess", "bfgs"};
38 
39 static const char *NTR_INIT[64] = {  "constant", "direction", "interpolation"};
40 
41 static const char *NTR_UPDATE[64] = {  "reduction", "interpolation"};
42 
43 /*  Routine for BFGS preconditioner */
44 static PetscErrorCode MatLMVMSolveShell(PC pc, Vec xin, Vec xout);
45 
46 /*
47    TaoSolve_NTR - Implements Newton's Method with a trust region approach
48    for solving unconstrained minimization problems.
49 
50    The basic algorithm is taken from MINPACK-2 (dstrn).
51 
52    TaoSolve_NTR computes a local minimizer of a twice differentiable function
53    f by applying a trust region variant of Newton's method.  At each stage
54    of the algorithm, we use the prconditioned conjugate gradient method to
55    determine an approximate minimizer of the quadratic equation
56 
57         q(s) = <s, Hs + g>
58 
59    subject to the trust region constraint
60 
61         || s ||_M <= radius,
62 
63    where radius is the trust region radius and M is a symmetric positive
64    definite matrix (the preconditioner).  Here g is the gradient and H
65    is the Hessian matrix.
66 
67    Note:  TaoSolve_NTR MUST use the iterative solver KSPNASH, KSPSTCG,
68           or KSPGLTR.  Thus, we set KSPNASH, KSPSTCG, or KSPGLTR in this
69           routine regardless of what the user may have previously specified.
70 */
71 #undef __FUNCT__
72 #define __FUNCT__ "TaoSolve_NTR"
73 static PetscErrorCode TaoSolve_NTR(Tao tao)
74 {
75   TAO_NTR            *tr = (TAO_NTR *)tao->data;
76   PC                 pc;
77   KSPConvergedReason ksp_reason;
78   TaoConvergedReason reason;
79   PetscReal          fmin, ftrial, prered, actred, kappa, sigma, beta;
80   PetscReal          tau, tau_1, tau_2, tau_max, tau_min, max_radius;
81   PetscReal          f, gnorm;
82 
83   PetscReal          delta;
84   PetscReal          norm_d;
85   PetscErrorCode     ierr;
86   PetscInt           iter = 0;
87   PetscInt           bfgsUpdates = 0;
88   PetscInt           needH;
89 
90   PetscInt           i_max = 5;
91   PetscInt           j_max = 1;
92   PetscInt           i, j, N, n, its;
93 
94   PetscFunctionBegin;
95   if (tao->XL || tao->XU || tao->ops->computebounds) {
96     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntr algorithm\n");CHKERRQ(ierr);
97   }
98 
99   tao->trust = tao->trust0;
100 
101   /* Modify the radius if it is too large or small */
102   tao->trust = PetscMax(tao->trust, tr->min_radius);
103   tao->trust = PetscMin(tao->trust, tr->max_radius);
104 
105 
106   if (NTR_PC_BFGS == tr->pc_type && !tr->M) {
107     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
108     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
109     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&tr->M);CHKERRQ(ierr);
110     ierr = MatLMVMAllocateVectors(tr->M,tao->solution);CHKERRQ(ierr);
111   }
112 
113   /* Check convergence criteria */
114   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
115   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
116   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
117   needH = 1;
118 
119   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
120   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
121 
122   /* Create vectors for the limited memory preconditioner */
123   if ((NTR_PC_BFGS == tr->pc_type) &&
124       (BFGS_SCALE_BFGS != tr->bfgs_scale_type)) {
125     if (!tr->Diag) {
126         ierr = VecDuplicate(tao->solution, &tr->Diag);CHKERRQ(ierr);
127     }
128   }
129 
130   switch(tr->ksp_type) {
131   case NTR_KSP_NASH:
132     ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr);
133     if (tao->ksp->ops->setfromoptions) {
134       (*tao->ksp->ops->setfromoptions)(tao->ksp);
135     }
136     break;
137 
138   case NTR_KSP_STCG:
139     ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr);
140     if (tao->ksp->ops->setfromoptions) {
141       (*tao->ksp->ops->setfromoptions)(tao->ksp);
142     }
143     break;
144 
145   default:
146     ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr);
147     if (tao->ksp->ops->setfromoptions) {
148       (*tao->ksp->ops->setfromoptions)(tao->ksp);
149     }
150     break;
151   }
152 
153   /*  Modify the preconditioner to use the bfgs approximation */
154   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
155   switch(tr->pc_type) {
156   case NTR_PC_NONE:
157     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
158     if (pc->ops->setfromoptions) {
159       (*pc->ops->setfromoptions)(pc);
160     }
161     break;
162 
163   case NTR_PC_AHESS:
164     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
165     if (pc->ops->setfromoptions) {
166       (*pc->ops->setfromoptions)(pc);
167     }
168     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
169     break;
170 
171   case NTR_PC_BFGS:
172     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
173     if (pc->ops->setfromoptions) {
174       (*pc->ops->setfromoptions)(pc);
175     }
176     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
177     ierr = PCShellSetContext(pc, tr->M);CHKERRQ(ierr);
178     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
179     break;
180 
181   default:
182     /*  Use the pc method set by pc_type */
183     break;
184   }
185 
186   /*  Initialize trust-region radius */
187   switch(tr->init_type) {
188   case NTR_INIT_CONSTANT:
189     /*  Use the initial radius specified */
190     break;
191 
192   case NTR_INIT_INTERPOLATION:
193     /*  Use the initial radius specified */
194     max_radius = 0.0;
195 
196     for (j = 0; j < j_max; ++j) {
197       fmin = f;
198       sigma = 0.0;
199 
200       if (needH) {
201         ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
202         needH = 0;
203       }
204 
205       for (i = 0; i < i_max; ++i) {
206 
207         ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
208         ierr = VecAXPY(tr->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr);
209         ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
210 
211         if (PetscIsInfOrNanReal(ftrial)) {
212           tau = tr->gamma1_i;
213         }
214         else {
215           if (ftrial < fmin) {
216             fmin = ftrial;
217             sigma = -tao->trust / gnorm;
218           }
219 
220           ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
221           ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
222 
223           prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm));
224           actred = f - ftrial;
225           if ((PetscAbsScalar(actred) <= tr->epsilon) &&
226               (PetscAbsScalar(prered) <= tr->epsilon)) {
227             kappa = 1.0;
228           }
229           else {
230             kappa = actred / prered;
231           }
232 
233           tau_1 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust + (1.0 - tr->theta_i) * prered - actred);
234           tau_2 = tr->theta_i * gnorm * tao->trust / (tr->theta_i * gnorm * tao->trust - (1.0 + tr->theta_i) * prered + actred);
235           tau_min = PetscMin(tau_1, tau_2);
236           tau_max = PetscMax(tau_1, tau_2);
237 
238           if (PetscAbsScalar(kappa - 1.0) <= tr->mu1_i) {
239             /*  Great agreement */
240             max_radius = PetscMax(max_radius, tao->trust);
241 
242             if (tau_max < 1.0) {
243               tau = tr->gamma3_i;
244             }
245             else if (tau_max > tr->gamma4_i) {
246               tau = tr->gamma4_i;
247             }
248             else {
249               tau = tau_max;
250             }
251           }
252           else if (PetscAbsScalar(kappa - 1.0) <= tr->mu2_i) {
253             /*  Good agreement */
254             max_radius = PetscMax(max_radius, tao->trust);
255 
256             if (tau_max < tr->gamma2_i) {
257               tau = tr->gamma2_i;
258             }
259             else if (tau_max > tr->gamma3_i) {
260               tau = tr->gamma3_i;
261             }
262             else {
263               tau = tau_max;
264             }
265           }
266           else {
267             /*  Not good agreement */
268             if (tau_min > 1.0) {
269               tau = tr->gamma2_i;
270             }
271             else if (tau_max < tr->gamma1_i) {
272               tau = tr->gamma1_i;
273             }
274             else if ((tau_min < tr->gamma1_i) && (tau_max >= 1.0)) {
275               tau = tr->gamma1_i;
276             }
277             else if ((tau_1 >= tr->gamma1_i) && (tau_1 < 1.0) &&
278                      ((tau_2 < tr->gamma1_i) || (tau_2 >= 1.0))) {
279               tau = tau_1;
280             }
281             else if ((tau_2 >= tr->gamma1_i) && (tau_2 < 1.0) &&
282                      ((tau_1 < tr->gamma1_i) || (tau_2 >= 1.0))) {
283               tau = tau_2;
284             }
285             else {
286               tau = tau_max;
287             }
288           }
289         }
290         tao->trust = tau * tao->trust;
291       }
292 
293       if (fmin < f) {
294         f = fmin;
295         ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr);
296         ierr = TaoComputeGradient(tao,tao->solution, tao->gradient);CHKERRQ(ierr);
297 
298         ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
299 
300         if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
301         needH = 1;
302 
303         ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr);
304         if (reason != TAO_CONTINUE_ITERATING) {
305           PetscFunctionReturn(0);
306         }
307       }
308     }
309     tao->trust = PetscMax(tao->trust, max_radius);
310 
311     /*  Modify the radius if it is too large or small */
312     tao->trust = PetscMax(tao->trust, tr->min_radius);
313     tao->trust = PetscMin(tao->trust, tr->max_radius);
314     break;
315 
316   default:
317     /*  Norm of the first direction will initialize radius */
318     tao->trust = 0.0;
319     break;
320   }
321 
322   /* Set initial scaling for the BFGS preconditioner
323      This step is done after computing the initial trust-region radius
324      since the function value may have decreased */
325   if (NTR_PC_BFGS == tr->pc_type) {
326     if (f != 0.0) {
327       delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
328     }
329     else {
330       delta = 2.0 / (gnorm*gnorm);
331     }
332     ierr = MatLMVMSetDelta(tr->M,delta);CHKERRQ(ierr);
333   }
334 
335   /* Have not converged; continue with Newton method */
336   while (reason == TAO_CONTINUE_ITERATING) {
337     ++iter;
338     tao->ksp_its=0;
339     /* Compute the Hessian */
340     if (needH) {
341       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
342       needH = 0;
343     }
344 
345     if (NTR_PC_BFGS == tr->pc_type) {
346       if (BFGS_SCALE_AHESS == tr->bfgs_scale_type) {
347         /* Obtain diagonal for the bfgs preconditioner */
348         ierr = MatGetDiagonal(tao->hessian, tr->Diag);CHKERRQ(ierr);
349         ierr = VecAbs(tr->Diag);CHKERRQ(ierr);
350         ierr = VecReciprocal(tr->Diag);CHKERRQ(ierr);
351         ierr = MatLMVMSetScale(tr->M,tr->Diag);CHKERRQ(ierr);
352       }
353 
354       /* Update the limited memory preconditioner */
355       ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
356       ++bfgsUpdates;
357     }
358 
359     while (reason == TAO_CONTINUE_ITERATING) {
360       ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr);
361 
362       /* Solve the trust region subproblem */
363       if (NTR_KSP_NASH == tr->ksp_type) {
364         ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
365         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
366         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
367         tao->ksp_its+=its;
368         tao->ksp_tot_its+=its;
369         ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
370       } else if (NTR_KSP_STCG == tr->ksp_type) {
371         ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
372         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
373         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
374         tao->ksp_its+=its;
375         tao->ksp_tot_its+=its;
376         ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
377       } else { /* NTR_KSP_GLTR */
378         ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
379         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
380         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
381         tao->ksp_its+=its;
382         tao->ksp_tot_its+=its;
383         ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
384       }
385 
386       if (0.0 == tao->trust) {
387         /* Radius was uninitialized; use the norm of the direction */
388         if (norm_d > 0.0) {
389           tao->trust = norm_d;
390 
391           /* Modify the radius if it is too large or small */
392           tao->trust = PetscMax(tao->trust, tr->min_radius);
393           tao->trust = PetscMin(tao->trust, tr->max_radius);
394         }
395         else {
396           /* The direction was bad; set radius to default value and re-solve
397              the trust-region subproblem to get a direction */
398           tao->trust = tao->trust0;
399 
400           /* Modify the radius if it is too large or small */
401           tao->trust = PetscMax(tao->trust, tr->min_radius);
402           tao->trust = PetscMin(tao->trust, tr->max_radius);
403 
404           if (NTR_KSP_NASH == tr->ksp_type) {
405             ierr = KSPNASHSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
406             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
407             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
408             tao->ksp_its+=its;
409             tao->ksp_tot_its+=its;
410             ierr = KSPNASHGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
411           } else if (NTR_KSP_STCG == tr->ksp_type) {
412             ierr = KSPSTCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
413             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
414             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
415             tao->ksp_its+=its;
416             tao->ksp_tot_its+=its;
417             ierr = KSPSTCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
418           } else { /* NTR_KSP_GLTR */
419             ierr = KSPGLTRSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
420             ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
421             ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
422             tao->ksp_its+=its;
423             tao->ksp_tot_its+=its;
424             ierr = KSPGLTRGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr);
425           }
426 
427           if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
428         }
429       }
430       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
431       ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
432       if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&
433           (NTR_PC_BFGS == tr->pc_type) && (bfgsUpdates > 1)) {
434         /* Preconditioner is numerically indefinite; reset the
435            approximate if using BFGS preconditioning. */
436 
437         if (f != 0.0) {
438           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
439         }
440         else {
441           delta = 2.0 / (gnorm*gnorm);
442         }
443         ierr = MatLMVMSetDelta(tr->M, delta);CHKERRQ(ierr);
444         ierr = MatLMVMReset(tr->M);CHKERRQ(ierr);
445         ierr = MatLMVMUpdate(tr->M, tao->solution, tao->gradient);CHKERRQ(ierr);
446         bfgsUpdates = 1;
447       }
448 
449       if (NTR_UPDATE_REDUCTION == tr->update_type) {
450         /* Get predicted reduction */
451         if (NTR_KSP_NASH == tr->ksp_type) {
452           ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
453         } else if (NTR_KSP_STCG == tr->ksp_type) {
454           ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
455         } else { /* gltr */
456           ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
457         }
458 
459         if (prered >= 0.0) {
460           /* The predicted reduction has the wrong sign.  This cannot
461              happen in infinite precision arithmetic.  Step should
462              be rejected! */
463           tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
464         }
465         else {
466           /* Compute trial step and function value */
467           ierr = VecCopy(tao->solution,tr->W);CHKERRQ(ierr);
468           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
469           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
470 
471           if (PetscIsInfOrNanReal(ftrial)) {
472             tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
473           } else {
474             /* Compute and actual reduction */
475             actred = f - ftrial;
476             prered = -prered;
477             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
478                 (PetscAbsScalar(prered) <= tr->epsilon)) {
479               kappa = 1.0;
480             }
481             else {
482               kappa = actred / prered;
483             }
484 
485             /* Accept or reject the step and update radius */
486             if (kappa < tr->eta1) {
487               /* Reject the step */
488               tao->trust = tr->alpha1 * PetscMin(tao->trust, norm_d);
489             }
490             else {
491               /* Accept the step */
492               if (kappa < tr->eta2) {
493                 /* Marginal bad step */
494                 tao->trust = tr->alpha2 * PetscMin(tao->trust, norm_d);
495               }
496               else if (kappa < tr->eta3) {
497                 /* Reasonable step */
498                 tao->trust = tr->alpha3 * tao->trust;
499               }
500               else if (kappa < tr->eta4) {
501                 /* Good step */
502                 tao->trust = PetscMax(tr->alpha4 * norm_d, tao->trust);
503               }
504               else {
505                 /* Very good step */
506                 tao->trust = PetscMax(tr->alpha5 * norm_d, tao->trust);
507               }
508               break;
509             }
510           }
511         }
512       }
513       else {
514         /* Get predicted reduction */
515         if (NTR_KSP_NASH == tr->ksp_type) {
516           ierr = KSPNASHGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
517         } else if (NTR_KSP_STCG == tr->ksp_type) {
518           ierr = KSPSTCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
519         } else { /* gltr */
520           ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
521         }
522 
523         if (prered >= 0.0) {
524           /* The predicted reduction has the wrong sign.  This cannot
525              happen in infinite precision arithmetic.  Step should
526              be rejected! */
527           tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
528         }
529         else {
530           ierr = VecCopy(tao->solution, tr->W);CHKERRQ(ierr);
531           ierr = VecAXPY(tr->W, 1.0, tao->stepdirection);CHKERRQ(ierr);
532           ierr = TaoComputeObjective(tao, tr->W, &ftrial);CHKERRQ(ierr);
533           if (PetscIsInfOrNanReal(ftrial)) {
534             tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
535           }
536           else {
537             ierr = VecDot(tao->gradient, tao->stepdirection, &beta);CHKERRQ(ierr);
538             actred = f - ftrial;
539             prered = -prered;
540             if ((PetscAbsScalar(actred) <= tr->epsilon) &&
541                 (PetscAbsScalar(prered) <= tr->epsilon)) {
542               kappa = 1.0;
543             }
544             else {
545               kappa = actred / prered;
546             }
547 
548             tau_1 = tr->theta * beta / (tr->theta * beta - (1.0 - tr->theta) * prered + actred);
549             tau_2 = tr->theta * beta / (tr->theta * beta + (1.0 + tr->theta) * prered - actred);
550             tau_min = PetscMin(tau_1, tau_2);
551             tau_max = PetscMax(tau_1, tau_2);
552 
553             if (kappa >= 1.0 - tr->mu1) {
554               /* Great agreement; accept step and update radius */
555               if (tau_max < 1.0) {
556                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
557               }
558               else if (tau_max > tr->gamma4) {
559                 tao->trust = PetscMax(tao->trust, tr->gamma4 * norm_d);
560               }
561               else {
562                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
563               }
564               break;
565             }
566             else if (kappa >= 1.0 - tr->mu2) {
567               /* Good agreement */
568 
569               if (tau_max < tr->gamma2) {
570                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
571               }
572               else if (tau_max > tr->gamma3) {
573                 tao->trust = PetscMax(tao->trust, tr->gamma3 * norm_d);
574               }
575               else if (tau_max < 1.0) {
576                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
577               }
578               else {
579                 tao->trust = PetscMax(tao->trust, tau_max * norm_d);
580               }
581               break;
582             }
583             else {
584               /* Not good agreement */
585               if (tau_min > 1.0) {
586                 tao->trust = tr->gamma2 * PetscMin(tao->trust, norm_d);
587               }
588               else if (tau_max < tr->gamma1) {
589                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
590               }
591               else if ((tau_min < tr->gamma1) && (tau_max >= 1.0)) {
592                 tao->trust = tr->gamma1 * PetscMin(tao->trust, norm_d);
593               }
594               else if ((tau_1 >= tr->gamma1) && (tau_1 < 1.0) &&
595                        ((tau_2 < tr->gamma1) || (tau_2 >= 1.0))) {
596                 tao->trust = tau_1 * PetscMin(tao->trust, norm_d);
597               }
598               else if ((tau_2 >= tr->gamma1) && (tau_2 < 1.0) &&
599                        ((tau_1 < tr->gamma1) || (tau_2 >= 1.0))) {
600                 tao->trust = tau_2 * PetscMin(tao->trust, norm_d);
601               }
602               else {
603                 tao->trust = tau_max * PetscMin(tao->trust, norm_d);
604               }
605             }
606           }
607         }
608       }
609 
610       /* The step computed was not good and the radius was decreased.
611          Monitor the radius to terminate. */
612       ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
613     }
614 
615     /* The radius may have been increased; modify if it is too large */
616     tao->trust = PetscMin(tao->trust, tr->max_radius);
617 
618     if (reason == TAO_CONTINUE_ITERATING) {
619       ierr = VecCopy(tr->W, tao->solution);CHKERRQ(ierr);
620       f = ftrial;
621       ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);
622       ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
623       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
624       needH = 1;
625       ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, tao->trust, &reason);CHKERRQ(ierr);
626     }
627   }
628   PetscFunctionReturn(0);
629 }
630 
631 /*------------------------------------------------------------*/
632 #undef __FUNCT__
633 #define __FUNCT__ "TaoSetUp_NTR"
634 static PetscErrorCode TaoSetUp_NTR(Tao tao)
635 {
636   TAO_NTR *tr = (TAO_NTR *)tao->data;
637   PetscErrorCode ierr;
638 
639   PetscFunctionBegin;
640 
641   if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);}
642   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);}
643   if (!tr->W) {ierr = VecDuplicate(tao->solution, &tr->W);CHKERRQ(ierr);}
644 
645   tr->Diag = 0;
646   tr->M = 0;
647 
648 
649   PetscFunctionReturn(0);
650 }
651 
652 /*------------------------------------------------------------*/
653 #undef __FUNCT__
654 #define __FUNCT__ "TaoDestroy_NTR"
655 static PetscErrorCode TaoDestroy_NTR(Tao tao)
656 {
657   TAO_NTR        *tr = (TAO_NTR *)tao->data;
658   PetscErrorCode ierr;
659 
660   PetscFunctionBegin;
661   if (tao->setupcalled) {
662     ierr = VecDestroy(&tr->W);CHKERRQ(ierr);
663   }
664   ierr = MatDestroy(&tr->M);CHKERRQ(ierr);
665   ierr = VecDestroy(&tr->Diag);CHKERRQ(ierr);
666   ierr = PetscFree(tao->data);CHKERRQ(ierr);
667   PetscFunctionReturn(0);
668 }
669 
670 /*------------------------------------------------------------*/
671 #undef __FUNCT__
672 #define __FUNCT__ "TaoSetFromOptions_NTR"
673 static PetscErrorCode TaoSetFromOptions_NTR(Tao tao)
674 {
675   TAO_NTR        *tr = (TAO_NTR *)tao->data;
676   PetscErrorCode ierr;
677 
678   PetscFunctionBegin;
679   ierr = PetscOptionsHead("Newton trust region method for unconstrained optimization");CHKERRQ(ierr);
680   ierr = PetscOptionsEList("-tao_ntr_ksp_type", "ksp type", "", NTR_KSP, NTR_KSP_TYPES, NTR_KSP[tr->ksp_type], &tr->ksp_type,NULL);CHKERRQ(ierr);
681   ierr = PetscOptionsEList("-tao_ntr_pc_type", "pc type", "", NTR_PC, NTR_PC_TYPES, NTR_PC[tr->pc_type], &tr->pc_type,NULL);CHKERRQ(ierr);
682   ierr = PetscOptionsEList("-tao_ntr_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[tr->bfgs_scale_type], &tr->bfgs_scale_type,NULL);CHKERRQ(ierr);
683   ierr = PetscOptionsEList("-tao_ntr_init_type", "tao->trust initialization type", "", NTR_INIT, NTR_INIT_TYPES, NTR_INIT[tr->init_type], &tr->init_type,NULL);CHKERRQ(ierr);
684   ierr = PetscOptionsEList("-tao_ntr_update_type", "radius update type", "", NTR_UPDATE, NTR_UPDATE_TYPES, NTR_UPDATE[tr->update_type], &tr->update_type,NULL);CHKERRQ(ierr);
685   ierr = PetscOptionsReal("-tao_ntr_eta1", "step is unsuccessful if actual reduction < eta1 * predicted reduction", "", tr->eta1, &tr->eta1,NULL);CHKERRQ(ierr);
686   ierr = PetscOptionsReal("-tao_ntr_eta2", "", "", tr->eta2, &tr->eta2,NULL);CHKERRQ(ierr);
687   ierr = PetscOptionsReal("-tao_ntr_eta3", "", "", tr->eta3, &tr->eta3,NULL);CHKERRQ(ierr);
688   ierr = PetscOptionsReal("-tao_ntr_eta4", "", "", tr->eta4, &tr->eta4,NULL);CHKERRQ(ierr);
689   ierr = PetscOptionsReal("-tao_ntr_alpha1", "", "", tr->alpha1, &tr->alpha1,NULL);CHKERRQ(ierr);
690   ierr = PetscOptionsReal("-tao_ntr_alpha2", "", "", tr->alpha2, &tr->alpha2,NULL);CHKERRQ(ierr);
691   ierr = PetscOptionsReal("-tao_ntr_alpha3", "", "", tr->alpha3, &tr->alpha3,NULL);CHKERRQ(ierr);
692   ierr = PetscOptionsReal("-tao_ntr_alpha4", "", "", tr->alpha4, &tr->alpha4,NULL);CHKERRQ(ierr);
693   ierr = PetscOptionsReal("-tao_ntr_alpha5", "", "", tr->alpha5, &tr->alpha5,NULL);CHKERRQ(ierr);
694   ierr = PetscOptionsReal("-tao_ntr_mu1", "", "", tr->mu1, &tr->mu1,NULL);CHKERRQ(ierr);
695   ierr = PetscOptionsReal("-tao_ntr_mu2", "", "", tr->mu2, &tr->mu2,NULL);CHKERRQ(ierr);
696   ierr = PetscOptionsReal("-tao_ntr_gamma1", "", "", tr->gamma1, &tr->gamma1,NULL);CHKERRQ(ierr);
697   ierr = PetscOptionsReal("-tao_ntr_gamma2", "", "", tr->gamma2, &tr->gamma2,NULL);CHKERRQ(ierr);
698   ierr = PetscOptionsReal("-tao_ntr_gamma3", "", "", tr->gamma3, &tr->gamma3,NULL);CHKERRQ(ierr);
699   ierr = PetscOptionsReal("-tao_ntr_gamma4", "", "", tr->gamma4, &tr->gamma4,NULL);CHKERRQ(ierr);
700   ierr = PetscOptionsReal("-tao_ntr_theta", "", "", tr->theta, &tr->theta,NULL);CHKERRQ(ierr);
701   ierr = PetscOptionsReal("-tao_ntr_mu1_i", "", "", tr->mu1_i, &tr->mu1_i,NULL);CHKERRQ(ierr);
702   ierr = PetscOptionsReal("-tao_ntr_mu2_i", "", "", tr->mu2_i, &tr->mu2_i,NULL);CHKERRQ(ierr);
703   ierr = PetscOptionsReal("-tao_ntr_gamma1_i", "", "", tr->gamma1_i, &tr->gamma1_i,NULL);CHKERRQ(ierr);
704   ierr = PetscOptionsReal("-tao_ntr_gamma2_i", "", "", tr->gamma2_i, &tr->gamma2_i,NULL);CHKERRQ(ierr);
705   ierr = PetscOptionsReal("-tao_ntr_gamma3_i", "", "", tr->gamma3_i, &tr->gamma3_i,NULL);CHKERRQ(ierr);
706   ierr = PetscOptionsReal("-tao_ntr_gamma4_i", "", "", tr->gamma4_i, &tr->gamma4_i,NULL);CHKERRQ(ierr);
707   ierr = PetscOptionsReal("-tao_ntr_theta_i", "", "", tr->theta_i, &tr->theta_i,NULL);CHKERRQ(ierr);
708   ierr = PetscOptionsReal("-tao_ntr_min_radius", "lower bound on initial trust-region radius", "", tr->min_radius, &tr->min_radius,NULL);CHKERRQ(ierr);
709   ierr = PetscOptionsReal("-tao_ntr_max_radius", "upper bound on trust-region radius", "", tr->max_radius, &tr->max_radius,NULL);CHKERRQ(ierr);
710   ierr = PetscOptionsReal("-tao_ntr_epsilon", "tolerance used when computing actual and predicted reduction", "", tr->epsilon, &tr->epsilon,NULL);CHKERRQ(ierr);
711   ierr = PetscOptionsTail();CHKERRQ(ierr);
712   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
713   PetscFunctionReturn(0);
714 }
715 
716 /*------------------------------------------------------------*/
717 #undef __FUNCT__
718 #define __FUNCT__ "TaoView_NTR"
719 static PetscErrorCode TaoView_NTR(Tao tao, PetscViewer viewer)
720 {
721   TAO_NTR        *tr = (TAO_NTR *)tao->data;
722   PetscErrorCode ierr;
723   PetscInt       nrejects;
724   PetscBool      isascii;
725 
726   PetscFunctionBegin;
727   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
728   if (isascii) {
729     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
730     if (NTR_PC_BFGS == tr->pc_type && tr->M) {
731       ierr = MatLMVMGetRejects(tr->M, &nrejects);CHKERRQ(ierr);
732       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n", nrejects);CHKERRQ(ierr);
733     }
734     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
735   }
736   PetscFunctionReturn(0);
737 }
738 
739 /*------------------------------------------------------------*/
740 /*MC
741   TAONTR - Newton's method with trust region for unconstrained minimization.
742   At each iteration, the Newton trust region method solves the system.
743   NTR expects a KSP solver with a trust region radius.
744             min_d  .5 dT Hk d + gkT d,  s.t.   ||d|| < Delta_k
745 
746   Options Database Keys:
747 + -tao_ntr_ksp_type - "nash","stcg","gltr"
748 . -tao_ntr_pc_type - "none","ahess","bfgs","petsc"
749 . -tao_ntr_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs"
750 . -tao_ntr_init_type - "constant","direction","interpolation"
751 . -tao_ntr_update_type - "reduction","interpolation"
752 . -tao_ntr_min_radius - lower bound on trust region radius
753 . -tao_ntr_max_radius - upper bound on trust region radius
754 . -tao_ntr_epsilon - tolerance for accepting actual / predicted reduction
755 . -tao_ntr_mu1_i - mu1 interpolation init factor
756 . -tao_ntr_mu2_i - mu2 interpolation init factor
757 . -tao_ntr_gamma1_i - gamma1 interpolation init factor
758 . -tao_ntr_gamma2_i - gamma2 interpolation init factor
759 . -tao_ntr_gamma3_i - gamma3 interpolation init factor
760 . -tao_ntr_gamma4_i - gamma4 interpolation init factor
761 . -tao_ntr_theta_i - thetha1 interpolation init factor
762 . -tao_ntr_eta1 - eta1 reduction update factor
763 . -tao_ntr_eta2 - eta2 reduction update factor
764 . -tao_ntr_eta3 - eta3 reduction update factor
765 . -tao_ntr_eta4 - eta4 reduction update factor
766 . -tao_ntr_alpha1 - alpha1 reduction update factor
767 . -tao_ntr_alpha2 - alpha2 reduction update factor
768 . -tao_ntr_alpha3 - alpha3 reduction update factor
769 . -tao_ntr_alpha4 - alpha4 reduction update factor
770 . -tao_ntr_alpha4 - alpha4 reduction update factor
771 . -tao_ntr_mu1 - mu1 interpolation update
772 . -tao_ntr_mu2 - mu2 interpolation update
773 . -tao_ntr_gamma1 - gamma1 interpolcation update
774 . -tao_ntr_gamma2 - gamma2 interpolcation update
775 . -tao_ntr_gamma3 - gamma3 interpolcation update
776 . -tao_ntr_gamma4 - gamma4 interpolation update
777 - -tao_ntr_theta - theta interpolation update
778 
779   Level: beginner
780 M*/
781 
782 #undef __FUNCT__
783 #define __FUNCT__ "TaoCreate_NTR"
784 PETSC_EXTERN PetscErrorCode TaoCreate_NTR(Tao tao)
785 {
786   TAO_NTR *tr;
787   PetscErrorCode ierr;
788 
789   PetscFunctionBegin;
790 
791   ierr = PetscNewLog(tao,&tr);CHKERRQ(ierr);
792 
793   tao->ops->setup = TaoSetUp_NTR;
794   tao->ops->solve = TaoSolve_NTR;
795   tao->ops->view = TaoView_NTR;
796   tao->ops->setfromoptions = TaoSetFromOptions_NTR;
797   tao->ops->destroy = TaoDestroy_NTR;
798 
799   tao->max_it = 50;
800 #if defined(PETSC_USE_REAL_SINGLE)
801   tao->fatol = 1e-5;
802   tao->frtol = 1e-5;
803 #else
804   tao->fatol = 1e-10;
805   tao->frtol = 1e-10;
806 #endif
807   tao->data = (void*)tr;
808 
809   tao->trust0 = 100.0;
810 
811   /*  Standard trust region update parameters */
812   tr->eta1 = 1.0e-4;
813   tr->eta2 = 0.25;
814   tr->eta3 = 0.50;
815   tr->eta4 = 0.90;
816 
817   tr->alpha1 = 0.25;
818   tr->alpha2 = 0.50;
819   tr->alpha3 = 1.00;
820   tr->alpha4 = 2.00;
821   tr->alpha5 = 4.00;
822 
823   /*  Interpolation parameters */
824   tr->mu1_i = 0.35;
825   tr->mu2_i = 0.50;
826 
827   tr->gamma1_i = 0.0625;
828   tr->gamma2_i = 0.50;
829   tr->gamma3_i = 2.00;
830   tr->gamma4_i = 5.00;
831 
832   tr->theta_i = 0.25;
833 
834   /*  Interpolation trust region update parameters */
835   tr->mu1 = 0.10;
836   tr->mu2 = 0.50;
837 
838   tr->gamma1 = 0.25;
839   tr->gamma2 = 0.50;
840   tr->gamma3 = 2.00;
841   tr->gamma4 = 4.00;
842 
843   tr->theta = 0.05;
844 
845   tr->min_radius = 1.0e-10;
846   tr->max_radius = 1.0e10;
847   tr->epsilon = 1.0e-6;
848 
849   tr->ksp_type        = NTR_KSP_STCG;
850   tr->pc_type         = NTR_PC_BFGS;
851   tr->bfgs_scale_type = BFGS_SCALE_AHESS;
852   tr->init_type       = NTR_INIT_INTERPOLATION;
853   tr->update_type     = NTR_UPDATE_REDUCTION;
854 
855 
856   /* Set linear solver to default for trust region */
857   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
858   PetscFunctionReturn(0);
859 }
860 
861 
862 #undef __FUNCT__
863 #define __FUNCT__ "MatLMVMSolveShell"
864 static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
865 {
866     PetscErrorCode ierr;
867     Mat M;
868     PetscFunctionBegin;
869     PetscValidHeaderSpecific(pc,PC_CLASSID,1);
870     PetscValidHeaderSpecific(b,VEC_CLASSID,2);
871     PetscValidHeaderSpecific(x,VEC_CLASSID,3);
872     ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
873     ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
874     PetscFunctionReturn(0);
875 }
876