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