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