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