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