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