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