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