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