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