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