xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 937a31a162199036c7d7383da36013b6a359c186)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/bound/impls/bnk/bnk.h>
3 
4 #include <petscksp.h>
5 
6 static const char *BNK_PC[64] = {"none", "ahess", "bfgs", "petsc"};
7 static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
8 static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
9 static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
10 static const char *BNK_AS[64] = {"none", "bertsekas"};
11 
12 /*------------------------------------------------------------*/
13 
14 /* Routine for BFGS preconditioner */
15 
16 PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
17 {
18   PetscErrorCode ierr;
19   Mat            M;
20 
21   PetscFunctionBegin;
22   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
23   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
24   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
25   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
26   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
27   PetscFunctionReturn(0);
28 }
29 
30 /*------------------------------------------------------------*/
31 
32 /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
33 
34 PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
35 {
36   PetscErrorCode               ierr;
37   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
38   PC                           pc;
39 
40   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma;
41   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
42   PetscReal                    resnorm, delta;
43 
44   PetscInt                     n,N;
45 
46   PetscInt                     i_max = 5;
47   PetscInt                     j_max = 1;
48   PetscInt                     i, j;
49 
50   PetscFunctionBegin;
51   /* Project the current point onto the feasible set */
52   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
53   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
54   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
55 
56   /* Project the initial point onto the feasible region */
57   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
58 
59   /* Check convergence criteria */
60   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
61   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
62   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
63   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
64 
65   /* Test the initial point for convergence */
66   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
67   ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
68   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
69   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
70   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
71   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
72 
73   /* Reset KSP stopping reason counters */
74   bnk->ksp_atol = 0;
75   bnk->ksp_rtol = 0;
76   bnk->ksp_dtol = 0;
77   bnk->ksp_ctol = 0;
78   bnk->ksp_negc = 0;
79   bnk->ksp_iter = 0;
80   bnk->ksp_othr = 0;
81 
82   /* Reset accepted step type counters */
83   bnk->tot_cg_its = 0;
84   bnk->newt = 0;
85   bnk->bfgs = 0;
86   bnk->sgrad = 0;
87   bnk->grad = 0;
88 
89   /* Initialize the Hessian perturbation */
90   bnk->pert = bnk->sval;
91 
92   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
93   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
94 
95   /* Allocate the vectors needed for the BFGS approximation */
96   if (BNK_PC_BFGS == bnk->pc_type) {
97     if (!bnk->M) {
98       ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
99       ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
100       ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
101       ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
102     }
103     if (bnk->bfgs_scale_type != BFGS_SCALE_BFGS && !bnk->Diag) {
104       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
105     }
106   }
107 
108   /* Prepare the min/max vectors for safeguarding diagonal scales */
109   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
110   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
111 
112   /* Modify the preconditioner to use the bfgs approximation */
113   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
114   switch(bnk->pc_type) {
115   case BNK_PC_NONE:
116     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
117     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
118     break;
119 
120   case BNK_PC_AHESS:
121     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
122     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
123     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
124     break;
125 
126   case BNK_PC_BFGS:
127     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
128     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
129     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
130     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
131     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
132     break;
133 
134   default:
135     /* Use the pc method set by pc_type */
136     break;
137   }
138 
139   /* Initialize trust-region radius.  The initialization is only performed
140      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
141   *needH = PETSC_TRUE;
142   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
143     switch(initType) {
144     case BNK_INIT_CONSTANT:
145       /* Use the initial radius specified */
146       tao->trust = tao->trust0;
147       break;
148 
149     case BNK_INIT_INTERPOLATION:
150       /* Use interpolation based on the initial Hessian */
151       max_radius = 0.0;
152       tao->trust = tao->trust0;
153       for (j = 0; j < j_max; ++j) {
154         f_min = bnk->f;
155         sigma = 0.0;
156 
157         if (*needH) {
158           /* Compute the Hessian at the new step, and extract the inactive subsystem */
159           ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
160           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
161           if (bnk->inactive_idx) {
162             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
163           } else {
164             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
165           }
166           *needH = PETSC_FALSE;
167         }
168 
169         for (i = 0; i < i_max; ++i) {
170           /* Take a steepest descent step and snap it to bounds */
171           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
172           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
173           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
174           /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */
175           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
176           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
177           /* Compute the objective at the trial */
178           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
179           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
180           if (PetscIsInfOrNanReal(ftrial)) {
181             tau = bnk->gamma1_i;
182           } else {
183             if (ftrial < f_min) {
184               f_min = ftrial;
185               sigma = -tao->trust / bnk->gnorm;
186             }
187 
188             /* Compute the predicted and actual reduction */
189             if (bnk->inactive_idx) {
190               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
191               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
192             } else {
193               bnk->X_inactive = bnk->W;
194               bnk->inactive_work = bnk->Xwork;
195             }
196             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
197             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
198             if (bnk->inactive_idx) {
199               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
200               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
201             }
202             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
203             actred = bnk->f - ftrial;
204             if ((PetscAbsScalar(actred) <= bnk->epsilon) &&
205                 (PetscAbsScalar(prered) <= bnk->epsilon)) {
206               kappa = 1.0;
207             }
208             else {
209               kappa = actred / prered;
210             }
211 
212             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
213             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
214             tau_min = PetscMin(tau_1, tau_2);
215             tau_max = PetscMax(tau_1, tau_2);
216 
217             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
218               /*  Great agreement */
219               max_radius = PetscMax(max_radius, tao->trust);
220 
221               if (tau_max < 1.0) {
222                 tau = bnk->gamma3_i;
223               }
224               else if (tau_max > bnk->gamma4_i) {
225                 tau = bnk->gamma4_i;
226               }
227               else {
228                 tau = tau_max;
229               }
230             }
231             else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
232               /*  Good agreement */
233               max_radius = PetscMax(max_radius, tao->trust);
234 
235               if (tau_max < bnk->gamma2_i) {
236                 tau = bnk->gamma2_i;
237               } else if (tau_max > bnk->gamma3_i) {
238                 tau = bnk->gamma3_i;
239               } else {
240                 tau = tau_max;
241               }
242             }
243             else {
244               /*  Not good agreement */
245               if (tau_min > 1.0) {
246                 tau = bnk->gamma2_i;
247               } else if (tau_max < bnk->gamma1_i) {
248                 tau = bnk->gamma1_i;
249               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
250                 tau = bnk->gamma1_i;
251               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) &&
252                         ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
253                 tau = tau_1;
254               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) &&
255                         ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
256                 tau = tau_2;
257               } else {
258                 tau = tau_max;
259               }
260             }
261           }
262           tao->trust = tau * tao->trust;
263         }
264 
265         if (f_min < bnk->f) {
266           /* We accidentally found a solution better than the initial, so accept it */
267           bnk->f = f_min;
268           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
269           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
270           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
271           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
272           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
273           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
274           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
275           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
276           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
277           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
278           *needH = PETSC_TRUE;
279           /* Test the new step for convergence */
280           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
281           ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
282           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
283           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
284           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
285           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
286           /* active BNCG recycling early because we have a stepdirection computed */
287           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
288         }
289       }
290       tao->trust = PetscMax(tao->trust, max_radius);
291 
292       /* Ensure that the trust radius is within the limits */
293       tao->trust = PetscMax(tao->trust, bnk->min_radius);
294       tao->trust = PetscMin(tao->trust, bnk->max_radius);
295       break;
296 
297     default:
298       /* Norm of the first direction will initialize radius */
299       tao->trust = 0.0;
300       break;
301     }
302   }
303 
304   /* Set initial scaling for the BFGS preconditioner
305      This step is done after computing the initial trust-region radius
306      since the function value may have decreased */
307   if (BNK_PC_BFGS == bnk->pc_type) {
308     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
309     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
310   }
311   PetscFunctionReturn(0);
312 }
313 
314 /*------------------------------------------------------------*/
315 
316 /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
317 
318 PetscErrorCode TaoBNKComputeHessian(Tao tao)
319 {
320   PetscErrorCode               ierr;
321   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
322 
323   PetscFunctionBegin;
324   /* Compute the Hessian */
325   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
326   /* Add a correction to the BFGS preconditioner */
327   if (BNK_PC_BFGS == bnk->pc_type) {
328     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
329     /* Update the BFGS diagonal scaling */
330     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type) {
331       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
332       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
333       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
334       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
335       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
336     }
337   }
338   PetscFunctionReturn(0);
339 }
340 
341 /*------------------------------------------------------------*/
342 
343 /* Routine for estimating the active set */
344 
345 PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
346 {
347   PetscErrorCode               ierr;
348   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
349 
350   PetscFunctionBegin;
351   switch (asType) {
352   case BNK_AS_NONE:
353     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
354     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
355     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
356     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
357     break;
358 
359   case BNK_AS_BERTSEKAS:
360     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
361     if (BNK_PC_BFGS == bnk->pc_type) {
362       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
363       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
364       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
365     } else {
366       /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
367       ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
368       ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
369       ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
370       ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
371       ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
372     }
373     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
374     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
375 
376   default:
377     break;
378   }
379   PetscFunctionReturn(0);
380 }
381 
382 /*------------------------------------------------------------*/
383 
384 /* Routine for bounding the step direction */
385 
386 PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step)
387 {
388   PetscErrorCode               ierr;
389   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
390 
391   PetscFunctionBegin;
392   if (bnk->active_idx) {
393     switch (bnk->as_type) {
394     case BNK_AS_NONE:
395       if (bnk->active_idx) {
396         ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
397         ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr);
398         ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
399       }
400       break;
401 
402     case BNK_AS_BERTSEKAS:
403       ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr);
404       break;
405 
406     default:
407       break;
408     }
409   }
410   PetscFunctionReturn(0);
411 }
412 
413 /*------------------------------------------------------------*/
414 
415 /* Routine for taking a finite number of BNCG iterations to
416    accelerate Newton convergence.
417 
418    In practice, this approach simply trades off Hessian evaluations
419    for more gradient evaluations.
420 */
421 
422 PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
423 {
424   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
425   PetscErrorCode               ierr;
426 
427   PetscFunctionBegin;
428   *terminate = PETSC_FALSE;
429   if (bnk->max_cg_its > 0) {
430     /* Copy the current solution, unprojected gradient and step info into BNCG */
431     bnk->bncg_ctx->f = bnk->f;
432     ierr = VecCopy(tao->solution, bnk->bncg->solution);CHKERRQ(ierr);
433     ierr = VecCopy(bnk->unprojected_gradient, bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
434     ierr = VecCopy(tao->stepdirection, bnk->bncg->stepdirection);CHKERRQ(ierr);
435     /* Take some small finite number of BNCG iterations */
436     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
437     /* Add the number of gradient and function evaluations to the total */
438     tao->nfuncs += bnk->bncg->nfuncs;
439     tao->nfuncgrads += bnk->bncg->nfuncgrads;
440     tao->ngrads += bnk->bncg->ngrads;
441     tao->nhess += bnk->bncg->nhess;
442     bnk->tot_cg_its += bnk->bncg->niter;
443     /* Extract the BNCG solution out and save it into BNK */
444     bnk->f = bnk->bncg_ctx->f;
445     ierr = VecCopy(bnk->bncg->solution, tao->solution);
446     ierr = VecCopy(bnk->bncg_ctx->unprojected_gradient, bnk->unprojected_gradient);
447     ierr = VecCopy(bnk->bncg->gradient, tao->gradient);CHKERRQ(ierr);
448     if (bnk->bncg->reason == TAO_CONVERGED_GATOL || bnk->bncg->reason == TAO_CONVERGED_GRTOL || bnk->bncg->reason == TAO_CONVERGED_GTTOL || bnk->bncg->reason == TAO_CONVERGED_MINF) {
449       *terminate = PETSC_TRUE;
450     }
451   }
452   PetscFunctionReturn(0);
453 }
454 
455 /*------------------------------------------------------------*/
456 
457 /* Routine for computing the Newton step. */
458 
459 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
460 {
461   PetscErrorCode               ierr;
462   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
463 
464   PetscReal                    delta;
465   PetscInt                     bfgsUpdates = 0;
466   PetscInt                     kspits;
467 
468   PetscFunctionBegin;
469   /* Prepare the reduced sub-matrices for the inactive set */
470   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
471   if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); }
472   if (bnk->inactive_idx) {
473     ierr = MatDestroy(&bnk->H_inactive);
474     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
475     if (tao->hessian == tao->hessian_pre) {
476       bnk->Hpre_inactive = bnk->H_inactive;
477     } else {
478       ierr = MatDestroy(&bnk->Hpre_inactive);
479       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
480     }
481   } else {
482     ierr = MatDestroy(&bnk->H_inactive);
483     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
484     if (tao->hessian == tao->hessian_pre) {
485       bnk->Hpre_inactive = bnk->H_inactive;
486     } else {
487       ierr = MatDestroy(&bnk->Hpre_inactive);
488       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
489     }
490   }
491 
492   /* Shift the reduced Hessian matrix */
493   if ((shift) && (bnk->pert > 0)) {
494     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
495     if (bnk->H_inactive != bnk->Hpre_inactive) {
496       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
497     }
498   }
499 
500   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
501   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
502     /* Obtain diagonal for the bfgs preconditioner  */
503     ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr);
504     if (bnk->inactive_idx) {
505       ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
506     } else {
507       bnk->Diag_red = bnk->Diag;
508     }
509     ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr);
510     if (bnk->inactive_idx) {
511       ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
512     }
513     ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
514     ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
515     ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
516     ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
517   }
518 
519   /* Solve the Newton system of equations */
520   tao->ksp_its = 0;
521   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
522   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
523   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
524   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
525   if (bnk->inactive_idx) {
526     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
527     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
528   } else {
529     bnk->G_inactive = bnk->unprojected_gradient;
530     bnk->X_inactive = tao->stepdirection;
531   }
532   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
533     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
534     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
535     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
536     tao->ksp_its+=kspits;
537     tao->ksp_tot_its+=kspits;
538     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
539 
540     if (0.0 == tao->trust) {
541       /* Radius was uninitialized; use the norm of the direction */
542       if (bnk->dnorm > 0.0) {
543         tao->trust = bnk->dnorm;
544 
545         /* Modify the radius if it is too large or small */
546         tao->trust = PetscMax(tao->trust, bnk->min_radius);
547         tao->trust = PetscMin(tao->trust, bnk->max_radius);
548       } else {
549         /* The direction was bad; set radius to default value and re-solve
550            the trust-region subproblem to get a direction */
551         tao->trust = tao->trust0;
552 
553         /* Modify the radius if it is too large or small */
554         tao->trust = PetscMax(tao->trust, bnk->min_radius);
555         tao->trust = PetscMin(tao->trust, bnk->max_radius);
556 
557         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
558         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
559         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
560         tao->ksp_its+=kspits;
561         tao->ksp_tot_its+=kspits;
562         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
563 
564         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
565       }
566     }
567   } else {
568     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
569     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
570     tao->ksp_its += kspits;
571     tao->ksp_tot_its+=kspits;
572   }
573   /* Restore sub vectors back */
574   if (bnk->inactive_idx) {
575     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
576     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
577   }
578   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
579   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
580   ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
581 
582   /* Record convergence reasons */
583   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
584   if (KSP_CONVERGED_ATOL == *ksp_reason) {
585     ++bnk->ksp_atol;
586   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
587     ++bnk->ksp_rtol;
588   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
589     ++bnk->ksp_ctol;
590   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
591     ++bnk->ksp_negc;
592   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
593     ++bnk->ksp_dtol;
594   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
595     ++bnk->ksp_iter;
596   } else {
597     ++bnk->ksp_othr;
598   }
599 
600   /* Make sure the BFGS preconditioner is healthy */
601   if (bnk->pc_type == BNK_PC_BFGS) {
602     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
603     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
604       /* Preconditioner is numerically indefinite; reset the approximation. */
605       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
606       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
607       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
608       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
609     }
610   }
611   PetscFunctionReturn(0);
612 }
613 
614 /*------------------------------------------------------------*/
615 
616 /* Routine for recomputing the predicted reduction for a given step vector */
617 
618 PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
619 {
620   PetscErrorCode               ierr;
621   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
622 
623   PetscFunctionBegin;
624   /* Extract subvectors associated with the inactive set */
625   if (bnk->inactive_idx){
626     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
627     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
628     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
629   } else {
630     bnk->X_inactive = tao->stepdirection;
631     bnk->inactive_work = bnk->Xwork;
632     bnk->G_inactive = bnk->Gwork;
633   }
634   /* Recompute the predicted decrease based on the quadratic model */
635   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
636   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
637   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);
638   /* Restore the sub vectors */
639   if (bnk->inactive_idx){
640     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
641     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
642     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
643   }
644   PetscFunctionReturn(0);
645 }
646 
647 /*------------------------------------------------------------*/
648 
649 /* Routine for ensuring that the Newton step is a descent direction.
650 
651    The step direction falls back onto BFGS, scaled gradient and gradient steps
652    in the event that the Newton step fails the test.
653 */
654 
655 PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
656 {
657   PetscErrorCode               ierr;
658   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
659 
660   PetscReal                    gdx, delta, e_min;
661   PetscInt                     bfgsUpdates;
662 
663   PetscFunctionBegin;
664   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
665   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
666     /* Newton step is not descent or direction produced Inf or NaN
667        Update the perturbation for next time */
668     if (bnk->pert <= 0.0) {
669       /* Initialize the perturbation */
670       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
671       if (bnk->is_gltr) {
672         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
673         bnk->pert = PetscMax(bnk->pert, -e_min);
674       }
675     } else {
676       /* Increase the perturbation */
677       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
678     }
679 
680     if (BNK_PC_BFGS != bnk->pc_type) {
681       /* We don't have the bfgs matrix around and updated
682          Must use gradient direction in this case */
683       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
684       *stepType = BNK_GRADIENT;
685     } else {
686       /* Attempt to use the BFGS direction */
687       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
688       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
689 
690       /* Check for success (descent direction)
691          NOTE: Negative gdx here means not a descent direction because
692          the fall-back step is missing a negative sign. */
693       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
694       if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
695         /* BFGS direction is not descent or direction produced not a number
696            We can assert bfgsUpdates > 1 in this case because
697            the first solve produces the scaled gradient direction,
698            which is guaranteed to be descent */
699 
700         /* Use steepest descent direction (scaled) */
701         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
702         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
703         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
704         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
705         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
706 
707         *stepType = BNK_SCALED_GRADIENT;
708       } else {
709         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
710         if (1 == bfgsUpdates) {
711           /* The first BFGS direction is always the scaled gradient */
712           *stepType = BNK_SCALED_GRADIENT;
713         } else {
714           *stepType = BNK_BFGS;
715         }
716       }
717     }
718     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
719     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
720     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
721   } else {
722     /* Computed Newton step is descent */
723     switch (ksp_reason) {
724     case KSP_DIVERGED_NANORINF:
725     case KSP_DIVERGED_BREAKDOWN:
726     case KSP_DIVERGED_INDEFINITE_MAT:
727     case KSP_DIVERGED_INDEFINITE_PC:
728     case KSP_CONVERGED_CG_NEG_CURVE:
729       /* Matrix or preconditioner is indefinite; increase perturbation */
730       if (bnk->pert <= 0.0) {
731         /* Initialize the perturbation */
732         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
733         if (bnk->is_gltr) {
734           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
735           bnk->pert = PetscMax(bnk->pert, -e_min);
736         }
737       } else {
738         /* Increase the perturbation */
739         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
740       }
741       break;
742 
743     default:
744       /* Newton step computation is good; decrease perturbation */
745       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
746       if (bnk->pert < bnk->pmin) {
747         bnk->pert = 0.0;
748       }
749       break;
750     }
751     *stepType = BNK_NEWTON;
752   }
753   PetscFunctionReturn(0);
754 }
755 
756 /*------------------------------------------------------------*/
757 
758 /* Routine for performing a bound-projected More-Thuente line search.
759 
760   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
761   Newton step does not produce a valid step length.
762 */
763 
764 PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
765 {
766   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
767   PetscErrorCode ierr;
768   TaoLineSearchConvergedReason ls_reason;
769 
770   PetscReal      e_min, gdx, delta;
771   PetscInt       bfgsUpdates;
772 
773   PetscFunctionBegin;
774   /* Perform the linesearch */
775   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
776   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
777 
778   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_GRADIENT) {
779     /* Linesearch failed, revert solution */
780     bnk->f = bnk->fold;
781     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
782     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
783 
784     switch(*stepType) {
785     case BNK_NEWTON:
786       /* Failed to obtain acceptable iterate with Newton step
787          Update the perturbation for next time */
788       if (bnk->pert <= 0.0) {
789         /* Initialize the perturbation */
790         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
791         if (bnk->is_gltr) {
792           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
793           bnk->pert = PetscMax(bnk->pert, -e_min);
794         }
795       } else {
796         /* Increase the perturbation */
797         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
798       }
799 
800       if (BNK_PC_BFGS != bnk->pc_type) {
801         /* We don't have the bfgs matrix around and being updated
802            Must use gradient direction in this case */
803         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
804         *stepType = BNK_GRADIENT;
805       } else {
806         /* Attempt to use the BFGS direction */
807         ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
808         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
809         /* Check for success (descent direction)
810            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
811         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
812         if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
813           /* BFGS direction is not descent or direction produced not a number
814              We can assert bfgsUpdates > 1 in this case
815              Use steepest descent direction (scaled) */
816           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
817           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
818           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
819           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
820           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
821 
822           bfgsUpdates = 1;
823           *stepType = BNK_SCALED_GRADIENT;
824         } else {
825           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
826           if (1 == bfgsUpdates) {
827             /* The first BFGS direction is always the scaled gradient */
828             *stepType = BNK_SCALED_GRADIENT;
829           } else {
830             *stepType = BNK_BFGS;
831           }
832         }
833       }
834       break;
835 
836     case BNK_BFGS:
837       /* Can only enter if pc_type == BNK_PC_BFGS
838          Failed to obtain acceptable iterate with BFGS step
839          Attempt to use the scaled gradient direction */
840       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
841       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
842       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
843       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
844       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
845       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
846 
847       bfgsUpdates = 1;
848       *stepType = BNK_SCALED_GRADIENT;
849       break;
850 
851     case BNK_SCALED_GRADIENT:
852       /* Can only enter if pc_type == BNK_PC_BFGS
853          The scaled gradient step did not produce a new iterate;
854          reset the BFGS matrix and attemp to use the gradient direction. */
855       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
856       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
857       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
858       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
859       ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
860 
861       bfgsUpdates = 1;
862       *stepType = BNK_GRADIENT;
863       break;
864     }
865     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
866     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
867     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
868 
869     /* Perform one last line search with the fall-back step */
870     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
871     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
872   }
873   *reason = ls_reason;
874   PetscFunctionReturn(0);
875 }
876 
877 /*------------------------------------------------------------*/
878 
879 /* Routine for updating the trust radius.
880 
881   Function features three different update methods:
882   1) Line-search step length based
883   2) Predicted decrease on the CG quadratic model
884   3) Interpolation
885 */
886 
887 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
888 {
889   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
890   PetscErrorCode ierr;
891 
892   PetscReal      step, kappa;
893   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
894 
895   PetscFunctionBegin;
896   /* Update trust region radius */
897   *accept = PETSC_FALSE;
898   switch(updateType) {
899   case BNK_UPDATE_STEP:
900     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
901     if (stepType == BNK_NEWTON) {
902       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
903       if (step < bnk->nu1) {
904         /* Very bad step taken; reduce radius */
905         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
906       } else if (step < bnk->nu2) {
907         /* Reasonably bad step taken; reduce radius */
908         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
909       } else if (step < bnk->nu3) {
910         /*  Reasonable step was taken; leave radius alone */
911         if (bnk->omega3 < 1.0) {
912           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
913         } else if (bnk->omega3 > 1.0) {
914           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
915         }
916       } else if (step < bnk->nu4) {
917         /*  Full step taken; increase the radius */
918         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
919       } else {
920         /*  More than full step taken; increase the radius */
921         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
922       }
923     } else {
924       /*  Newton step was not good; reduce the radius */
925       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
926     }
927     break;
928 
929   case BNK_UPDATE_REDUCTION:
930     if (stepType == BNK_NEWTON) {
931       if (prered < 0.0) {
932         /* The predicted reduction has the wrong sign.  This cannot
933            happen in infinite precision arithmetic.  Step should
934            be rejected! */
935         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
936       }
937       else {
938         if (PetscIsInfOrNanReal(actred)) {
939           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
940         } else {
941           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) &&
942               (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
943             kappa = 1.0;
944           }
945           else {
946             kappa = actred / prered;
947           }
948 
949           /* Accept or reject the step and update radius */
950           if (kappa < bnk->eta1) {
951             /* Reject the step */
952             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
953           }
954           else {
955             /* Accept the step */
956             *accept = PETSC_TRUE;
957             /* Update the trust region radius only if the computed step is at the trust radius boundary */
958             if (bnk->dnorm == tao->trust) {
959               if (kappa < bnk->eta2) {
960                 /* Marginal bad step */
961                 tao->trust = bnk->alpha2 * tao->trust;
962               }
963               else if (kappa < bnk->eta3) {
964                 /* Reasonable step */
965                 tao->trust = bnk->alpha3 * tao->trust;
966               }
967               else if (kappa < bnk->eta4) {
968                 /* Good step */
969                 tao->trust = bnk->alpha4 * tao->trust;
970               }
971               else {
972                 /* Very good step */
973                 tao->trust = bnk->alpha5 * tao->trust;
974               }
975             }
976           }
977         }
978       }
979     } else {
980       /*  Newton step was not good; reduce the radius */
981       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
982     }
983     break;
984 
985   default:
986     if (stepType == BNK_NEWTON) {
987       if (prered < 0.0) {
988         /*  The predicted reduction has the wrong sign.  This cannot */
989         /*  happen in infinite precision arithmetic.  Step should */
990         /*  be rejected! */
991         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
992       } else {
993         if (PetscIsInfOrNanReal(actred)) {
994           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
995         } else {
996           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
997             kappa = 1.0;
998           } else {
999             kappa = actred / prered;
1000           }
1001 
1002           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
1003           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
1004           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
1005           tau_min = PetscMin(tau_1, tau_2);
1006           tau_max = PetscMax(tau_1, tau_2);
1007 
1008           if (kappa >= 1.0 - bnk->mu1) {
1009             /*  Great agreement */
1010             *accept = PETSC_TRUE;
1011             if (tau_max < 1.0) {
1012               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1013             } else if (tau_max > bnk->gamma4) {
1014               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
1015             } else {
1016               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1017             }
1018           } else if (kappa >= 1.0 - bnk->mu2) {
1019             /*  Good agreement */
1020             *accept = PETSC_TRUE;
1021             if (tau_max < bnk->gamma2) {
1022               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1023             } else if (tau_max > bnk->gamma3) {
1024               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1025             } else if (tau_max < 1.0) {
1026               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1027             } else {
1028               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1029             }
1030           } else {
1031             /*  Not good agreement */
1032             if (tau_min > 1.0) {
1033               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1034             } else if (tau_max < bnk->gamma1) {
1035               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1036             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
1037               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1038             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
1039               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
1040             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
1041               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
1042             } else {
1043               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1044             }
1045           }
1046         }
1047       }
1048     } else {
1049       /*  Newton step was not good; reduce the radius */
1050       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
1051     }
1052     break;
1053   }
1054   /* Make sure the radius does not violate min and max settings */
1055   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1056   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1057   PetscFunctionReturn(0);
1058 }
1059 
1060 /* ---------------------------------------------------------- */
1061 
1062 PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
1063 {
1064   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1065 
1066   PetscFunctionBegin;
1067   switch (stepType) {
1068   case BNK_NEWTON:
1069     ++bnk->newt;
1070     break;
1071   case BNK_BFGS:
1072     ++bnk->bfgs;
1073     break;
1074   case BNK_SCALED_GRADIENT:
1075     ++bnk->sgrad;
1076     break;
1077   case BNK_GRADIENT:
1078     ++bnk->grad;
1079     break;
1080   default:
1081     break;
1082   }
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 /* ---------------------------------------------------------- */
1087 
1088 PetscErrorCode TaoSetUp_BNK(Tao tao)
1089 {
1090   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1091   PetscErrorCode ierr;
1092   KSPType        ksp_type;
1093   PetscInt       i;
1094 
1095   PetscFunctionBegin;
1096   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
1097   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
1098   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
1099   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
1100   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
1101   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
1102   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
1103   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
1104   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
1105   if (!bnk->Diag_min) {ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);}
1106   if (!bnk->Diag_max) {ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);}
1107   if (bnk->max_cg_its > 0) {
1108     if (!bnk->bncg_sol) {ierr = VecDuplicate(tao->solution,&bnk->bncg_sol);CHKERRQ(ierr);}
1109     ierr = TaoSetInitialVector(bnk->bncg, bnk->bncg_sol);CHKERRQ(ierr);
1110     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1111     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1112     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1113 
1114     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1115     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1116     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1117     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1118     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1119     for (i=0; i<tao->numbermonitors; i++) {
1120       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1121       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1122     }
1123     ierr = TaoSetUp(bnk->bncg);CHKERRQ(ierr);
1124     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
1125   }
1126   bnk->Diag = 0;
1127   bnk->Diag_red = 0;
1128   bnk->X_inactive = 0;
1129   bnk->G_inactive = 0;
1130   bnk->inactive_work = 0;
1131   bnk->active_work = 0;
1132   bnk->inactive_idx = 0;
1133   bnk->active_idx = 0;
1134   bnk->active_lower = 0;
1135   bnk->active_upper = 0;
1136   bnk->active_fixed = 0;
1137   bnk->M = 0;
1138   bnk->H_inactive = 0;
1139   bnk->Hpre_inactive = 0;
1140   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
1141   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
1142   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
1143   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1144   PetscFunctionReturn(0);
1145 }
1146 
1147 /*------------------------------------------------------------*/
1148 
1149 static PetscErrorCode TaoDestroy_BNK(Tao tao)
1150 {
1151   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1152   PetscErrorCode ierr;
1153 
1154   PetscFunctionBegin;
1155   if (tao->setupcalled) {
1156     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1157     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1158     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
1159     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
1160     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
1161     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
1162     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1163     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
1164     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1165     if (bnk->max_cg_its > 0) {
1166       ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1167       ierr = VecDestroy(&bnk->bncg_sol);CHKERRQ(ierr);
1168     }
1169   }
1170   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1171   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1172   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);}
1173   if (bnk->H_inactive != tao->hessian) {ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);}
1174   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 /*------------------------------------------------------------*/
1179 
1180 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1181 {
1182   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1183   PetscErrorCode ierr;
1184 
1185   PetscFunctionBegin;
1186   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
1187   ierr = PetscOptionsEList("-tao_bnk_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr);
1188   ierr = PetscOptionsEList("-tao_bnk_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[bnk->bfgs_scale_type], &bnk->bfgs_scale_type, 0);CHKERRQ(ierr);
1189   ierr = PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, 0);CHKERRQ(ierr);
1190   ierr = PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);CHKERRQ(ierr);
1191   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr);
1192   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1193   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1194   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1195   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1196   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1197   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1198   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1199   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1200   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1201   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1202   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1203   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1204   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1205   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1206   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
1207   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
1208   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
1209   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
1210   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
1211   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
1212   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
1213   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
1214   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
1215   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
1216   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
1217   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
1218   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
1219   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
1220   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
1221   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
1222   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
1223   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
1224   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
1225   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
1226   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
1227   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
1228   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1229   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
1230   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
1231   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
1232   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
1233   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1234   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1235   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1236   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1237   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
1238   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1239   ierr = PetscOptionsInt("-tao_bnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its,NULL);CHKERRQ(ierr);
1240   ierr = PetscOptionsTail();CHKERRQ(ierr);
1241   ierr = TaoSetFromOptions(bnk->bncg);
1242   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1243   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1244   PetscFunctionReturn(0);
1245 }
1246 
1247 /*------------------------------------------------------------*/
1248 
1249 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1250 {
1251   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1252   PetscInt       nrejects;
1253   PetscBool      isascii;
1254   PetscErrorCode ierr;
1255 
1256   PetscFunctionBegin;
1257   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1258   if (isascii) {
1259     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1260     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1261       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1262       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1263     }
1264     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1265     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1266     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1267     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1268     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1269     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1270     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1271     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1272     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1273     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1274     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1275     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1276     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1277     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1278   }
1279   PetscFunctionReturn(0);
1280 }
1281 
1282 /* ---------------------------------------------------------- */
1283 
1284 /*MC
1285   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
1286   At each iteration, the BNK methods solve the symmetric
1287   system of equations to obtain the step diretion dk:
1288               Hk dk = -gk
1289   for free variables only. The step can be globalized either through
1290   trust-region methods, or a line search, or a heuristic mixture of both.
1291 
1292     Options Database Keys:
1293 + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
1294 . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
1295 . -tao_bnk_init_type - "constant","direction","interpolation"
1296 . -tao_bnk_update_type - "step","direction","interpolation"
1297 . -tao_bnk_as_type - "none","bertsekas"
1298 . -tao_bnk_sval - perturbation starting value
1299 . -tao_bnk_imin - minimum initial perturbation
1300 . -tao_bnk_imax - maximum initial perturbation
1301 . -tao_bnk_pmin - minimum perturbation
1302 . -tao_bnk_pmax - maximum perturbation
1303 . -tao_bnk_pgfac - growth factor
1304 . -tao_bnk_psfac - shrink factor
1305 . -tao_bnk_imfac - initial merit factor
1306 . -tao_bnk_pmgfac - merit growth factor
1307 . -tao_bnk_pmsfac - merit shrink factor
1308 . -tao_bnk_eta1 - poor steplength; reduce radius
1309 . -tao_bnk_eta2 - reasonable steplength; leave radius
1310 . -tao_bnk_eta3 - good steplength; increase readius
1311 . -tao_bnk_eta4 - excellent steplength; greatly increase radius
1312 . -tao_bnk_alpha1 - alpha1 reduction
1313 . -tao_bnk_alpha2 - alpha2 reduction
1314 . -tao_bnk_alpha3 - alpha3 reduction
1315 . -tao_bnk_alpha4 - alpha4 reduction
1316 . -tao_bnk_alpha - alpha5 reduction
1317 . -tao_bnk_mu1 - mu1 interpolation update
1318 . -tao_bnk_mu2 - mu2 interpolation update
1319 . -tao_bnk_gamma1 - gamma1 interpolation update
1320 . -tao_bnk_gamma2 - gamma2 interpolation update
1321 . -tao_bnk_gamma3 - gamma3 interpolation update
1322 . -tao_bnk_gamma4 - gamma4 interpolation update
1323 . -tao_bnk_theta - theta interpolation update
1324 . -tao_bnk_omega1 - omega1 step update
1325 . -tao_bnk_omega2 - omega2 step update
1326 . -tao_bnk_omega3 - omega3 step update
1327 . -tao_bnk_omega4 - omega4 step update
1328 . -tao_bnk_omega5 - omega5 step update
1329 . -tao_bnk_mu1_i -  mu1 interpolation init factor
1330 . -tao_bnk_mu2_i -  mu2 interpolation init factor
1331 . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
1332 . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
1333 . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
1334 . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
1335 . -tao_bnk_theta_i -  theta interpolation init factor
1336 - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1337 
1338   Level: beginner
1339 M*/
1340 
1341 PetscErrorCode TaoCreate_BNK(Tao tao)
1342 {
1343   TAO_BNK        *bnk;
1344   const char     *morethuente_type = TAOLINESEARCHMT;
1345   PetscErrorCode ierr;
1346 
1347   PetscFunctionBegin;
1348   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1349 
1350   tao->ops->setup = TaoSetUp_BNK;
1351   tao->ops->view = TaoView_BNK;
1352   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1353   tao->ops->destroy = TaoDestroy_BNK;
1354 
1355   /*  Override default settings (unless already changed) */
1356   if (!tao->max_it_changed) tao->max_it = 50;
1357   if (!tao->trust0_changed) tao->trust0 = 100.0;
1358 
1359   tao->data = (void*)bnk;
1360 
1361   /*  Hessian shifting parameters */
1362   bnk->sval   = 0.0;
1363   bnk->imin   = 1.0e-4;
1364   bnk->imax   = 1.0e+2;
1365   bnk->imfac  = 1.0e-1;
1366 
1367   bnk->pmin   = 1.0e-12;
1368   bnk->pmax   = 1.0e+2;
1369   bnk->pgfac  = 1.0e+1;
1370   bnk->psfac  = 4.0e-1;
1371   bnk->pmgfac = 1.0e-1;
1372   bnk->pmsfac = 1.0e-1;
1373 
1374   /*  Default values for trust-region radius update based on steplength */
1375   bnk->nu1 = 0.25;
1376   bnk->nu2 = 0.50;
1377   bnk->nu3 = 1.00;
1378   bnk->nu4 = 1.25;
1379 
1380   bnk->omega1 = 0.25;
1381   bnk->omega2 = 0.50;
1382   bnk->omega3 = 1.00;
1383   bnk->omega4 = 2.00;
1384   bnk->omega5 = 4.00;
1385 
1386   /*  Default values for trust-region radius update based on reduction */
1387   bnk->eta1 = 1.0e-4;
1388   bnk->eta2 = 0.25;
1389   bnk->eta3 = 0.50;
1390   bnk->eta4 = 0.90;
1391 
1392   bnk->alpha1 = 0.25;
1393   bnk->alpha2 = 0.50;
1394   bnk->alpha3 = 1.00;
1395   bnk->alpha4 = 2.00;
1396   bnk->alpha5 = 4.00;
1397 
1398   /*  Default values for trust-region radius update based on interpolation */
1399   bnk->mu1 = 0.10;
1400   bnk->mu2 = 0.50;
1401 
1402   bnk->gamma1 = 0.25;
1403   bnk->gamma2 = 0.50;
1404   bnk->gamma3 = 2.00;
1405   bnk->gamma4 = 4.00;
1406 
1407   bnk->theta = 0.05;
1408 
1409   /*  Default values for trust region initialization based on interpolation */
1410   bnk->mu1_i = 0.35;
1411   bnk->mu2_i = 0.50;
1412 
1413   bnk->gamma1_i = 0.0625;
1414   bnk->gamma2_i = 0.5;
1415   bnk->gamma3_i = 2.0;
1416   bnk->gamma4_i = 5.0;
1417 
1418   bnk->theta_i = 0.25;
1419 
1420   /*  Remaining parameters */
1421   bnk->max_cg_its = 0;
1422   bnk->min_radius = 1.0e-10;
1423   bnk->max_radius = 1.0e10;
1424   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
1425   bnk->as_tol = 1.0e-3;
1426   bnk->as_step = 1.0e-3;
1427   bnk->dmin = 1.0e-6;
1428   bnk->dmax = 1.0e6;
1429 
1430   bnk->pc_type         = BNK_PC_AHESS;
1431   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1432   bnk->init_type       = BNK_INIT_INTERPOLATION;
1433   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
1434   bnk->as_type         = BNK_AS_BERTSEKAS;
1435 
1436   /* Create the embedded BNCG solver */
1437   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1438   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1439   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1440   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1441 
1442   /* Create the line search */
1443   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1444   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1445   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1446   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1447   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1448 
1449   /*  Set linear solver to default for symmetric matrices */
1450   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1451   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1452   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1453   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1454   PetscFunctionReturn(0);
1455 }
1456