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