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