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