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