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