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