xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision b2d8c57778fe5daf492d4a621e5fa4036873e062)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener #include <petscksp.h>
5eb910715SAlp Dener 
6e031d6f5SAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
7e031d6f5SAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
8e031d6f5SAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
9e031d6f5SAlp Dener 
10e031d6f5SAlp Dener /*------------------------------------------------------------*/
11e031d6f5SAlp Dener 
12cd929ea3SAlp Dener PetscErrorCode TaoBNKPreconBFGS(PC BFGSpc, Vec X, Vec Y)
13eb910715SAlp Dener {
14eb910715SAlp Dener   PetscErrorCode ierr;
15cd929ea3SAlp Dener   Mat *M;
16eb910715SAlp Dener 
17eb910715SAlp Dener   PetscFunctionBegin;
18cd929ea3SAlp Dener   ierr = PCShellGetContext(BFGSpc, (void**)&M);CHKERRQ(ierr);
19cd929ea3SAlp Dener   ierr = MatSolve(*M, X, Y);CHKERRQ(ierr);
20eb910715SAlp Dener   PetscFunctionReturn(0);
21eb910715SAlp Dener }
22eb910715SAlp Dener 
23df278d8fSAlp Dener /*------------------------------------------------------------*/
24df278d8fSAlp Dener 
25df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
26df278d8fSAlp Dener 
27c0f10754SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
28eb910715SAlp Dener {
29eb910715SAlp Dener   PetscErrorCode               ierr;
30eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
31eb910715SAlp Dener   PC                           pc;
32eb910715SAlp Dener 
3389da521bSAlp Dener   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma, resnorm;
34eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
35b9ac7092SAlp Dener   PetscBool                    is_bfgs, is_jacobi;
36c4b75bccSAlp Dener   PetscInt                     n, N, nDiff;
37eb910715SAlp Dener   PetscInt                     i_max = 5;
38eb910715SAlp Dener   PetscInt                     j_max = 1;
39eb910715SAlp Dener   PetscInt                     i, j;
40eb910715SAlp Dener 
41eb910715SAlp Dener   PetscFunctionBegin;
4228017e9fSAlp Dener   /* Project the current point onto the feasible set */
4328017e9fSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
44e031d6f5SAlp Dener   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
45b9ac7092SAlp Dener   if (tao->bounded) {
4628017e9fSAlp Dener     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
47cd929ea3SAlp Dener   }
4828017e9fSAlp Dener 
4928017e9fSAlp Dener   /* Project the initial point onto the feasible region */
503b063059SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
5128017e9fSAlp Dener 
5228017e9fSAlp Dener   /* Check convergence criteria */
5328017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
5461be54a6SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
5561be54a6SAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
5661be54a6SAlp Dener   ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
5708752603SAlp Dener   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
5828017e9fSAlp Dener 
59c0f10754SAlp Dener   /* Test the initial point for convergence */
6089da521bSAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
6189da521bSAlp Dener   ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
62b4a30f08SAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
63e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
64e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
65c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
66c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
67c0f10754SAlp Dener 
68e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
69eb910715SAlp Dener   bnk->ksp_atol = 0;
70eb910715SAlp Dener   bnk->ksp_rtol = 0;
71eb910715SAlp Dener   bnk->ksp_dtol = 0;
72eb910715SAlp Dener   bnk->ksp_ctol = 0;
73eb910715SAlp Dener   bnk->ksp_negc = 0;
74eb910715SAlp Dener   bnk->ksp_iter = 0;
75eb910715SAlp Dener   bnk->ksp_othr = 0;
76eb910715SAlp Dener 
77e031d6f5SAlp Dener   /* Reset accepted step type counters */
78e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
79e031d6f5SAlp Dener   bnk->newt = 0;
80e031d6f5SAlp Dener   bnk->bfgs = 0;
81e031d6f5SAlp Dener   bnk->sgrad = 0;
82e031d6f5SAlp Dener   bnk->grad = 0;
83e031d6f5SAlp Dener 
84fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
85fed79b8eSAlp Dener   bnk->pert = bnk->sval;
86fed79b8eSAlp Dener 
87937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
88937a31a1SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
89937a31a1SAlp Dener 
90e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
91b9ac7092SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
92b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr);
93b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr);
94b9ac7092SAlp Dener   if (is_bfgs) {
95b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
96b9ac7092SAlp Dener     ierr = PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M);CHKERRQ(ierr);
97eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
98eb910715SAlp Dener     ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
99b9ac7092SAlp Dener     ierr = MatSetSizes(bnk->M, n, n, N, N);CHKERRQ(ierr);
100cd929ea3SAlp Dener     ierr = MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
101b9ac7092SAlp Dener   } else if (is_jacobi) {
102b9ac7092SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
103e031d6f5SAlp Dener   }
104e031d6f5SAlp Dener 
105e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
10662675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
10762675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
108eb910715SAlp Dener 
109eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
110eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
111c0f10754SAlp Dener   *needH = PETSC_TRUE;
112eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
11362675beeSAlp Dener     switch(initType) {
114eb910715SAlp Dener     case BNK_INIT_CONSTANT:
115eb910715SAlp Dener       /* Use the initial radius specified */
116c0f10754SAlp Dener       tao->trust = tao->trust0;
117eb910715SAlp Dener       break;
118eb910715SAlp Dener 
119eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
120c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
121eb910715SAlp Dener       max_radius = 0.0;
12208752603SAlp Dener       tao->trust = tao->trust0;
123eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1240a4511e9SAlp Dener         f_min = bnk->f;
125eb910715SAlp Dener         sigma = 0.0;
126eb910715SAlp Dener 
127c0f10754SAlp Dener         if (*needH) {
12862602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
129e031d6f5SAlp Dener           ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
13008752603SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
13189da521bSAlp Dener           ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
13289da521bSAlp Dener           if (bnk->active_idx) {
1332ab2a32cSAlp Dener             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
13428017e9fSAlp Dener           } else {
13508752603SAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
13628017e9fSAlp Dener           }
137c0f10754SAlp Dener           *needH = PETSC_FALSE;
138eb910715SAlp Dener         }
139eb910715SAlp Dener 
140eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
14162602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
14262602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
14362602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
1443b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
14589da521bSAlp Dener           /* Compute the step we actually accepted */
146eb910715SAlp Dener           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
14762602cfbSAlp Dener           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
14862602cfbSAlp Dener           /* Compute the objective at the trial */
14962602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
150b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
15162602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
152eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
153eb910715SAlp Dener             tau = bnk->gamma1_i;
154eb910715SAlp Dener           } else {
1550a4511e9SAlp Dener             if (ftrial < f_min) {
1560a4511e9SAlp Dener               f_min = ftrial;
157eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
158eb910715SAlp Dener             }
15908752603SAlp Dener 
160770b7498SAlp Dener             /* Compute the predicted and actual reduction */
16189da521bSAlp Dener             if (bnk->active_idx) {
16208752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
16308752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1642ab2a32cSAlp Dener             } else {
16508752603SAlp Dener               bnk->X_inactive = bnk->W;
16608752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1672ab2a32cSAlp Dener             }
16808752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
16908752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
17089da521bSAlp Dener             if (bnk->active_idx) {
17108752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
17208752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1732ab2a32cSAlp Dener             }
174eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
175eb910715SAlp Dener             actred = bnk->f - ftrial;
1763105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
177eb910715SAlp Dener               kappa = 1.0;
1783105154fSTodd Munson             } else {
179eb910715SAlp Dener               kappa = actred / prered;
180eb910715SAlp Dener             }
181eb910715SAlp Dener 
182eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
183eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
184eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
185eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
186eb910715SAlp Dener 
187eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
188eb910715SAlp Dener               /*  Great agreement */
189eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
190eb910715SAlp Dener 
191eb910715SAlp Dener               if (tau_max < 1.0) {
192eb910715SAlp Dener                 tau = bnk->gamma3_i;
1933105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
194eb910715SAlp Dener                 tau = bnk->gamma4_i;
1953105154fSTodd Munson               } else {
196eb910715SAlp Dener                 tau = tau_max;
197eb910715SAlp Dener               }
1988f8a4e06SAlp Dener             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
199eb910715SAlp Dener               /*  Good agreement */
200eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
201eb910715SAlp Dener 
202eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
203eb910715SAlp Dener                 tau = bnk->gamma2_i;
204eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
205eb910715SAlp Dener                 tau = bnk->gamma3_i;
206eb910715SAlp Dener               } else {
207eb910715SAlp Dener                 tau = tau_max;
208eb910715SAlp Dener               }
2098f8a4e06SAlp Dener             } else {
210eb910715SAlp Dener               /*  Not good agreement */
211eb910715SAlp Dener               if (tau_min > 1.0) {
212eb910715SAlp Dener                 tau = bnk->gamma2_i;
213eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
214eb910715SAlp Dener                 tau = bnk->gamma1_i;
215eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
216eb910715SAlp Dener                 tau = bnk->gamma1_i;
2173105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
218eb910715SAlp Dener                 tau = tau_1;
2193105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
220eb910715SAlp Dener                 tau = tau_2;
221eb910715SAlp Dener               } else {
222eb910715SAlp Dener                 tau = tau_max;
223eb910715SAlp Dener               }
224eb910715SAlp Dener             }
225eb910715SAlp Dener           }
226eb910715SAlp Dener           tao->trust = tau * tao->trust;
227eb910715SAlp Dener         }
228eb910715SAlp Dener 
2290a4511e9SAlp Dener         if (f_min < bnk->f) {
230937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2310a4511e9SAlp Dener           bnk->f = f_min;
232937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
233eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
2343b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
235937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
236937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
23709164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
23861be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
23961be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
24061be54a6SAlp Dener           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
241937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
242c4b75bccSAlp Dener           ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
243c0f10754SAlp Dener           *needH = PETSC_TRUE;
244937a31a1SAlp Dener           /* Test the new step for convergence */
24589da521bSAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
24689da521bSAlp Dener           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
247b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
248e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
249e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
250eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
251eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
252937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
253937a31a1SAlp Dener           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
254eb910715SAlp Dener         }
255eb910715SAlp Dener       }
256eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
257e031d6f5SAlp Dener 
258e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
259e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
260e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
261eb910715SAlp Dener       break;
262eb910715SAlp Dener 
263eb910715SAlp Dener     default:
264eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
265eb910715SAlp Dener       tao->trust = 0.0;
266eb910715SAlp Dener       break;
267eb910715SAlp Dener     }
268eb910715SAlp Dener   }
269eb910715SAlp Dener   PetscFunctionReturn(0);
270eb910715SAlp Dener }
271eb910715SAlp Dener 
272df278d8fSAlp Dener /*------------------------------------------------------------*/
273df278d8fSAlp Dener 
27462675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
27562675beeSAlp Dener 
27662675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
27762675beeSAlp Dener {
27862675beeSAlp Dener   PetscErrorCode               ierr;
27962675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
28062675beeSAlp Dener 
28162675beeSAlp Dener   PetscFunctionBegin;
28262675beeSAlp Dener   /* Compute the Hessian */
28362675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
28462675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
285b9ac7092SAlp Dener   if (bnk->M) {
28662675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
28762675beeSAlp Dener   }
28862675beeSAlp Dener   PetscFunctionReturn(0);
28962675beeSAlp Dener }
29062675beeSAlp Dener 
29162675beeSAlp Dener /*------------------------------------------------------------*/
29262675beeSAlp Dener 
2932f75a4aaSAlp Dener /* Routine for estimating the active set */
2942f75a4aaSAlp Dener 
29508752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
2962f75a4aaSAlp Dener {
2972f75a4aaSAlp Dener   PetscErrorCode               ierr;
2982f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
299c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3002f75a4aaSAlp Dener 
3012f75a4aaSAlp Dener   PetscFunctionBegin;
30208752603SAlp Dener   switch (asType) {
3032f75a4aaSAlp Dener   case BNK_AS_NONE:
3042f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3052f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3062f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3072f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3082f75a4aaSAlp Dener     break;
3092f75a4aaSAlp Dener 
3102f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3112f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
312b9ac7092SAlp Dener     if (bnk->M) {
3132f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
314cd929ea3SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3152f75a4aaSAlp Dener     } else {
31661be54a6SAlp Dener       ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
317c4b75bccSAlp Dener       ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
318c4b75bccSAlp Dener       if (hessComputed && diagExists) {
3199b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3202f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3219b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3229b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3232f75a4aaSAlp Dener         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3242f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
32561be54a6SAlp Dener       } else {
326c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
32761be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
32861be54a6SAlp Dener       }
3292f75a4aaSAlp Dener     }
3302f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
33189da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
33289da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
333c4b75bccSAlp Dener     break;
3342f75a4aaSAlp Dener 
3352f75a4aaSAlp Dener   default:
3362f75a4aaSAlp Dener     break;
3372f75a4aaSAlp Dener   }
3382f75a4aaSAlp Dener   PetscFunctionReturn(0);
3392f75a4aaSAlp Dener }
3402f75a4aaSAlp Dener 
3412f75a4aaSAlp Dener /*------------------------------------------------------------*/
3422f75a4aaSAlp Dener 
3432f75a4aaSAlp Dener /* Routine for bounding the step direction */
3442f75a4aaSAlp Dener 
345a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
3462f75a4aaSAlp Dener {
3472f75a4aaSAlp Dener   PetscErrorCode               ierr;
3482f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3492f75a4aaSAlp Dener 
3502f75a4aaSAlp Dener   PetscFunctionBegin;
351a1318120SAlp Dener   switch (asType) {
3522f75a4aaSAlp Dener   case BNK_AS_NONE:
353c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
3542f75a4aaSAlp Dener     break;
3552f75a4aaSAlp Dener 
3562f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
357c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
3582f75a4aaSAlp Dener     break;
3592f75a4aaSAlp Dener 
3602f75a4aaSAlp Dener   default:
3612f75a4aaSAlp Dener     break;
3622f75a4aaSAlp Dener   }
3632f75a4aaSAlp Dener   PetscFunctionReturn(0);
3642f75a4aaSAlp Dener }
3652f75a4aaSAlp Dener 
366e031d6f5SAlp Dener /*------------------------------------------------------------*/
367e031d6f5SAlp Dener 
368e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
369e031d6f5SAlp Dener    accelerate Newton convergence.
370e031d6f5SAlp Dener 
371e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
372e031d6f5SAlp Dener    for more gradient evaluations.
373e031d6f5SAlp Dener */
374e031d6f5SAlp Dener 
375c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
376c0f10754SAlp Dener {
377c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
378c0f10754SAlp Dener   PetscErrorCode               ierr;
379c0f10754SAlp Dener 
380c0f10754SAlp Dener   PetscFunctionBegin;
381c0f10754SAlp Dener   *terminate = PETSC_FALSE;
382c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
383c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
384c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
385c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
386c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
387c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
388c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
389c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
390c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
391c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
392e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
393c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
394c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
395c0f10754SAlp Dener     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) {
396c0f10754SAlp Dener       *terminate = PETSC_TRUE;
39761be54a6SAlp Dener     } else {
39833c78596SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
399c0f10754SAlp Dener     }
400c0f10754SAlp Dener   }
401c0f10754SAlp Dener   PetscFunctionReturn(0);
402c0f10754SAlp Dener }
403c0f10754SAlp Dener 
4042f75a4aaSAlp Dener /*------------------------------------------------------------*/
4052f75a4aaSAlp Dener 
406c0f10754SAlp Dener /* Routine for computing the Newton step. */
407df278d8fSAlp Dener 
40862675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
409eb910715SAlp Dener {
410eb910715SAlp Dener   PetscErrorCode               ierr;
411eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
412eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
413eb910715SAlp Dener   PetscInt                     kspits;
414eb910715SAlp Dener 
415eb910715SAlp Dener   PetscFunctionBegin;
41689da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
41789da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
41889da521bSAlp Dener   if (!bnk->inactive_idx) {
41989da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
420a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
42189da521bSAlp Dener     PetscFunctionReturn(0);
42289da521bSAlp Dener   }
42389da521bSAlp Dener 
4245e9b73cbSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
42589da521bSAlp Dener   if (bnk->active_idx) {
42633c78596SAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
4275e9b73cbSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
428eb3ba6a7SAlp Dener     if (tao->hessian == tao->hessian_pre) {
429eb3ba6a7SAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
430eb3ba6a7SAlp Dener     } else {
43133c78596SAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
4325e9b73cbSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
433eb3ba6a7SAlp Dener     }
434b9ac7092SAlp Dener     if (bnk->bfgs_pre) {
435b9ac7092SAlp Dener       ierr = PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx);CHKERRQ(ierr);
436b9ac7092SAlp Dener     }
4372f75a4aaSAlp Dener   } else {
43833c78596SAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
43933c78596SAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
44062675beeSAlp Dener     if (tao->hessian == tao->hessian_pre) {
44162675beeSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
44262675beeSAlp Dener     } else {
44333c78596SAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
44433c78596SAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);CHKERRQ(ierr);
44562675beeSAlp Dener     }
446b9ac7092SAlp Dener     if (bnk->bfgs_pre) {
447b9ac7092SAlp Dener       ierr = PCLMVMClearIS(bnk->bfgs_pre);CHKERRQ(ierr);
448b9ac7092SAlp Dener     }
44962675beeSAlp Dener   }
45062675beeSAlp Dener 
45162675beeSAlp Dener   /* Shift the reduced Hessian matrix */
45262675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
45362675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
45462675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
45562675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
45662675beeSAlp Dener     }
45762675beeSAlp Dener   }
45862675beeSAlp Dener 
459eb910715SAlp Dener   /* Solve the Newton system of equations */
460937a31a1SAlp Dener   tao->ksp_its = 0;
4612f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
4625e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
46309164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
4645e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
46589da521bSAlp Dener   if (bnk->active_idx) {
4665e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
4675e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
4685e9b73cbSAlp Dener   } else {
4695e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4705e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
47128017e9fSAlp Dener   }
472eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
473fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4745e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
475eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
476eb910715SAlp Dener     tao->ksp_its+=kspits;
477eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
478080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
479eb910715SAlp Dener 
480eb910715SAlp Dener     if (0.0 == tao->trust) {
481eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
482080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
483080d2917SAlp Dener         tao->trust = bnk->dnorm;
484eb910715SAlp Dener 
485eb910715SAlp Dener         /* Modify the radius if it is too large or small */
486eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
487eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
488eb910715SAlp Dener       } else {
489eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
490eb910715SAlp Dener            the trust-region subproblem to get a direction */
491eb910715SAlp Dener         tao->trust = tao->trust0;
492eb910715SAlp Dener 
493eb910715SAlp Dener         /* Modify the radius if it is too large or small */
494eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
495eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
496eb910715SAlp Dener 
497fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4985e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
499eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
500eb910715SAlp Dener         tao->ksp_its+=kspits;
501eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
502080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
503eb910715SAlp Dener 
504080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
505eb910715SAlp Dener       }
506eb910715SAlp Dener     }
507eb910715SAlp Dener   } else {
5085e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
509eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
510eb910715SAlp Dener     tao->ksp_its += kspits;
511eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
512eb910715SAlp Dener   }
5135e9b73cbSAlp Dener   /* Restore sub vectors back */
51489da521bSAlp Dener   if (bnk->active_idx) {
5155e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5165e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5175e9b73cbSAlp Dener   }
518770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
519fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
520a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
521770b7498SAlp Dener 
522770b7498SAlp Dener   /* Record convergence reasons */
523e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
524e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
525770b7498SAlp Dener     ++bnk->ksp_atol;
526e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
527770b7498SAlp Dener     ++bnk->ksp_rtol;
528e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
529770b7498SAlp Dener     ++bnk->ksp_ctol;
530e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
531770b7498SAlp Dener     ++bnk->ksp_negc;
532e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
533770b7498SAlp Dener     ++bnk->ksp_dtol;
534e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
535770b7498SAlp Dener     ++bnk->ksp_iter;
536770b7498SAlp Dener   } else {
537770b7498SAlp Dener     ++bnk->ksp_othr;
538770b7498SAlp Dener   }
539fed79b8eSAlp Dener 
540fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
541b9ac7092SAlp Dener   if (bnk->M) {
542cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
543*b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
544fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
545*b2d8c577SAlp Dener       ierr = MatLMVMResetJ0(bnk->M);CHKERRQ(ierr);
546cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
54709164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
548eb910715SAlp Dener     }
549fed79b8eSAlp Dener   }
550e465cd6fSAlp Dener   PetscFunctionReturn(0);
551e465cd6fSAlp Dener }
552eb910715SAlp Dener 
55362675beeSAlp Dener /*------------------------------------------------------------*/
55462675beeSAlp Dener 
5555e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5565e9b73cbSAlp Dener 
5575e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
5585e9b73cbSAlp Dener {
5595e9b73cbSAlp Dener   PetscErrorCode               ierr;
5605e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
5615e9b73cbSAlp Dener 
5625e9b73cbSAlp Dener   PetscFunctionBegin;
5635e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
56489da521bSAlp Dener   if (bnk->active_idx){
5655e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5665e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5675e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5685e9b73cbSAlp Dener   } else {
5695e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
5705e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5715e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
5725e9b73cbSAlp Dener   }
5735e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5745e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
5755e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
57633c78596SAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);CHKERRQ(ierr);
5775e9b73cbSAlp Dener   /* Restore the sub vectors */
57889da521bSAlp Dener   if (bnk->active_idx){
5795e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5805e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5815e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5825e9b73cbSAlp Dener   }
5835e9b73cbSAlp Dener   PetscFunctionReturn(0);
5845e9b73cbSAlp Dener }
5855e9b73cbSAlp Dener 
5865e9b73cbSAlp Dener /*------------------------------------------------------------*/
5875e9b73cbSAlp Dener 
58862675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
58962675beeSAlp Dener 
59062675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
59162675beeSAlp Dener    in the event that the Newton step fails the test.
59262675beeSAlp Dener */
59362675beeSAlp Dener 
594e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
595e465cd6fSAlp Dener {
596e465cd6fSAlp Dener   PetscErrorCode               ierr;
597e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
598e465cd6fSAlp Dener 
599*b2d8c577SAlp Dener   PetscReal                    gdx, e_min;
600e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
601e465cd6fSAlp Dener 
602e465cd6fSAlp Dener   PetscFunctionBegin;
603080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
604eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
605eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
606eb910715SAlp Dener        Update the perturbation for next time */
607eb910715SAlp Dener     if (bnk->pert <= 0.0) {
608eb910715SAlp Dener       /* Initialize the perturbation */
609eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
610eb910715SAlp Dener       if (bnk->is_gltr) {
611eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
612eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
613eb910715SAlp Dener       }
614eb910715SAlp Dener     } else {
615eb910715SAlp Dener       /* Increase the perturbation */
616eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
617eb910715SAlp Dener     }
618eb910715SAlp Dener 
619b9ac7092SAlp Dener     if (bnk->M) {
620eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
621eb910715SAlp Dener          Must use gradient direction in this case */
622080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
623eb910715SAlp Dener       *stepType = BNK_GRADIENT;
624eb910715SAlp Dener     } else {
625eb910715SAlp Dener       /* Attempt to use the BFGS direction */
626cd929ea3SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
627eb910715SAlp Dener 
6288d5ead36SAlp Dener       /* Check for success (descent direction)
6298d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
6308d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
631080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
6323105154fSTodd Munson       if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
633eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
634eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
635eb910715SAlp Dener            the first solve produces the scaled gradient direction,
636eb910715SAlp Dener            which is guaranteed to be descent */
637eb910715SAlp Dener 
638eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
639*b2d8c577SAlp Dener         ierr = MatLMVMResetJ0(bnk->M);CHKERRQ(ierr);
640cd929ea3SAlp Dener         ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
64109164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
642cd929ea3SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
643eb910715SAlp Dener 
644eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
645eb910715SAlp Dener       } else {
646cd929ea3SAlp Dener         ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
647eb910715SAlp Dener         if (1 == bfgsUpdates) {
648eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
649eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
650eb910715SAlp Dener         } else {
651eb910715SAlp Dener           *stepType = BNK_BFGS;
652eb910715SAlp Dener         }
653eb910715SAlp Dener       }
654eb910715SAlp Dener     }
6558d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6568d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
657a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
658eb910715SAlp Dener   } else {
659eb910715SAlp Dener     /* Computed Newton step is descent */
660eb910715SAlp Dener     switch (ksp_reason) {
661eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
662eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
663eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
664eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
665eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
666eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
667eb910715SAlp Dener       if (bnk->pert <= 0.0) {
668eb910715SAlp Dener         /* Initialize the perturbation */
669eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
670eb910715SAlp Dener         if (bnk->is_gltr) {
671eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
672eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
673eb910715SAlp Dener         }
674eb910715SAlp Dener       } else {
675eb910715SAlp Dener         /* Increase the perturbation */
676eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
677eb910715SAlp Dener       }
678eb910715SAlp Dener       break;
679eb910715SAlp Dener 
680eb910715SAlp Dener     default:
681eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
682eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
683eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
684eb910715SAlp Dener         bnk->pert = 0.0;
685eb910715SAlp Dener       }
686eb910715SAlp Dener       break;
687eb910715SAlp Dener     }
688fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
689eb910715SAlp Dener   }
690eb910715SAlp Dener   PetscFunctionReturn(0);
691eb910715SAlp Dener }
692eb910715SAlp Dener 
693df278d8fSAlp Dener /*------------------------------------------------------------*/
694df278d8fSAlp Dener 
695df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
696df278d8fSAlp Dener 
697df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
698df278d8fSAlp Dener   Newton step does not produce a valid step length.
699df278d8fSAlp Dener */
700df278d8fSAlp Dener 
701937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
702c14b763aSAlp Dener {
703c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
704c14b763aSAlp Dener   PetscErrorCode ierr;
705c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
706c14b763aSAlp Dener 
707*b2d8c577SAlp Dener   PetscReal      e_min, gdx;
708c14b763aSAlp Dener   PetscInt       bfgsUpdates;
709c14b763aSAlp Dener 
710c14b763aSAlp Dener   PetscFunctionBegin;
711c14b763aSAlp Dener   /* Perform the linesearch */
712c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
713c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
714c14b763aSAlp Dener 
715*b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
716c14b763aSAlp Dener     /* Linesearch failed, revert solution */
717c14b763aSAlp Dener     bnk->f = bnk->fold;
718c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
719c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
720c14b763aSAlp Dener 
721937a31a1SAlp Dener     switch(*stepType) {
722c14b763aSAlp Dener     case BNK_NEWTON:
7238d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
724c14b763aSAlp Dener          Update the perturbation for next time */
725c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
726c14b763aSAlp Dener         /* Initialize the perturbation */
727c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
728c14b763aSAlp Dener         if (bnk->is_gltr) {
729c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
730c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
731c14b763aSAlp Dener         }
732c14b763aSAlp Dener       } else {
733c14b763aSAlp Dener         /* Increase the perturbation */
734c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
735c14b763aSAlp Dener       }
736c14b763aSAlp Dener 
737b9ac7092SAlp Dener       if (bnk->M) {
738c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
739c14b763aSAlp Dener            Must use gradient direction in this case */
740937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
741937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
742c14b763aSAlp Dener       } else {
743c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
744cd929ea3SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
7458d5ead36SAlp Dener         /* Check for success (descent direction)
7468d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
747c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7483105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
749c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
750c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
751c14b763aSAlp Dener              Use steepest descent direction (scaled) */
752*b2d8c577SAlp Dener           ierr = MatLMVMResetJ0(bnk->M);CHKERRQ(ierr);
753cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
754c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
755cd929ea3SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
756c14b763aSAlp Dener 
757c14b763aSAlp Dener           bfgsUpdates = 1;
758937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
759c14b763aSAlp Dener         } else {
760cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
761c14b763aSAlp Dener           if (1 == bfgsUpdates) {
762c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
763937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
764c14b763aSAlp Dener           } else {
765937a31a1SAlp Dener             *stepType = BNK_BFGS;
766c14b763aSAlp Dener           }
767c14b763aSAlp Dener         }
768c14b763aSAlp Dener       }
769c14b763aSAlp Dener       break;
770c14b763aSAlp Dener 
771c14b763aSAlp Dener     case BNK_BFGS:
772c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
773c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
774c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
775*b2d8c577SAlp Dener       ierr = MatLMVMResetJ0(bnk->M);CHKERRQ(ierr);
776cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
777c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
778cd929ea3SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
779c14b763aSAlp Dener 
780c14b763aSAlp Dener       bfgsUpdates = 1;
781937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
782c14b763aSAlp Dener       break;
783c14b763aSAlp Dener     }
7848d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7858d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
786a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
787c14b763aSAlp Dener 
7888d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
789c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
790c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
791c14b763aSAlp Dener   }
792c14b763aSAlp Dener   *reason = ls_reason;
793c14b763aSAlp Dener   PetscFunctionReturn(0);
794c14b763aSAlp Dener }
795c14b763aSAlp Dener 
796df278d8fSAlp Dener /*------------------------------------------------------------*/
797df278d8fSAlp Dener 
798df278d8fSAlp Dener /* Routine for updating the trust radius.
799df278d8fSAlp Dener 
800df278d8fSAlp Dener   Function features three different update methods:
801df278d8fSAlp Dener   1) Line-search step length based
802df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
803df278d8fSAlp Dener   3) Interpolation
804df278d8fSAlp Dener */
805df278d8fSAlp Dener 
80628017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
807080d2917SAlp Dener {
808080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
809080d2917SAlp Dener   PetscErrorCode ierr;
810080d2917SAlp Dener 
811b1c2d0e3SAlp Dener   PetscReal      step, kappa;
812080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
813080d2917SAlp Dener 
814080d2917SAlp Dener   PetscFunctionBegin;
815080d2917SAlp Dener   /* Update trust region radius */
816080d2917SAlp Dener   *accept = PETSC_FALSE;
81728017e9fSAlp Dener   switch(updateType) {
818080d2917SAlp Dener   case BNK_UPDATE_STEP:
819c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
820080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
821080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
822080d2917SAlp Dener       if (step < bnk->nu1) {
823080d2917SAlp Dener         /* Very bad step taken; reduce radius */
824080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
825080d2917SAlp Dener       } else if (step < bnk->nu2) {
826080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
827080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
828080d2917SAlp Dener       } else if (step < bnk->nu3) {
829080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
830080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
831080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
832080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
833080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
834080d2917SAlp Dener         }
835080d2917SAlp Dener       } else if (step < bnk->nu4) {
836080d2917SAlp Dener         /*  Full step taken; increase the radius */
837080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
838080d2917SAlp Dener       } else {
839080d2917SAlp Dener         /*  More than full step taken; increase the radius */
840080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
841080d2917SAlp Dener       }
842080d2917SAlp Dener     } else {
843080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
844080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
845080d2917SAlp Dener     }
846080d2917SAlp Dener     break;
847080d2917SAlp Dener 
848080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
849080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
850b1c2d0e3SAlp Dener       if (prered < 0.0) {
851fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
852fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
853fed79b8eSAlp Dener            be rejected! */
854080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8553105154fSTodd Munson       } else {
856b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
857080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
858080d2917SAlp Dener         } else {
8593105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
860080d2917SAlp Dener             kappa = 1.0;
8613105154fSTodd Munson           } else {
862080d2917SAlp Dener             kappa = actred / prered;
863080d2917SAlp Dener           }
864fed79b8eSAlp Dener 
865fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
866080d2917SAlp Dener           if (kappa < bnk->eta1) {
867fed79b8eSAlp Dener             /* Reject the step */
868080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8693105154fSTodd Munson           } else {
870fed79b8eSAlp Dener             /* Accept the step */
871c133c014SAlp Dener             *accept = PETSC_TRUE;
872c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8738d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
874080d2917SAlp Dener               if (kappa < bnk->eta2) {
875080d2917SAlp Dener                 /* Marginal bad step */
876c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
8773105154fSTodd Munson               } else if (kappa < bnk->eta3) {
878fed79b8eSAlp Dener                 /* Reasonable step */
879fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
8803105154fSTodd Munson               } else if (kappa < bnk->eta4) {
881080d2917SAlp Dener                 /* Good step */
882c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
8833105154fSTodd Munson               } else {
884080d2917SAlp Dener                 /* Very good step */
885c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
886080d2917SAlp Dener               }
887c133c014SAlp Dener             }
888080d2917SAlp Dener           }
889080d2917SAlp Dener         }
890080d2917SAlp Dener       }
891080d2917SAlp Dener     } else {
892080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
893080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
894080d2917SAlp Dener     }
895080d2917SAlp Dener     break;
896080d2917SAlp Dener 
897080d2917SAlp Dener   default:
898080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
899b1c2d0e3SAlp Dener       if (prered < 0.0) {
900080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
901080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
902080d2917SAlp Dener         /*  be rejected! */
903080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
904080d2917SAlp Dener       } else {
905b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
906080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
907080d2917SAlp Dener         } else {
908080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
909080d2917SAlp Dener             kappa = 1.0;
910080d2917SAlp Dener           } else {
911080d2917SAlp Dener             kappa = actred / prered;
912080d2917SAlp Dener           }
913080d2917SAlp Dener 
914080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
915080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
916080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
917080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
918080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
919080d2917SAlp Dener 
920080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
921080d2917SAlp Dener             /*  Great agreement */
922080d2917SAlp Dener             *accept = PETSC_TRUE;
923080d2917SAlp Dener             if (tau_max < 1.0) {
924080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
925080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
926080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
927080d2917SAlp Dener             } else {
928080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
929080d2917SAlp Dener             }
930080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
931080d2917SAlp Dener             /*  Good agreement */
932080d2917SAlp Dener             *accept = PETSC_TRUE;
933080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
934080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
935080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
936080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
937080d2917SAlp Dener             } else if (tau_max < 1.0) {
938080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
939080d2917SAlp Dener             } else {
940080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
941080d2917SAlp Dener             }
942080d2917SAlp Dener           } else {
943080d2917SAlp Dener             /*  Not good agreement */
944080d2917SAlp Dener             if (tau_min > 1.0) {
945080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
946080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
947080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
948080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
949080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
950080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
951080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
952080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
953080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
954080d2917SAlp Dener             } else {
955080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
956080d2917SAlp Dener             }
957080d2917SAlp Dener           }
958080d2917SAlp Dener         }
959080d2917SAlp Dener       }
960080d2917SAlp Dener     } else {
961080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
962080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
963080d2917SAlp Dener     }
96428017e9fSAlp Dener     break;
965080d2917SAlp Dener   }
966c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
967c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
968fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
969080d2917SAlp Dener   PetscFunctionReturn(0);
970080d2917SAlp Dener }
971080d2917SAlp Dener 
972eb910715SAlp Dener /* ---------------------------------------------------------- */
973df278d8fSAlp Dener 
97462675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
97562675beeSAlp Dener {
97662675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
97762675beeSAlp Dener 
97862675beeSAlp Dener   PetscFunctionBegin;
97962675beeSAlp Dener   switch (stepType) {
98062675beeSAlp Dener   case BNK_NEWTON:
98162675beeSAlp Dener     ++bnk->newt;
98262675beeSAlp Dener     break;
98362675beeSAlp Dener   case BNK_BFGS:
98462675beeSAlp Dener     ++bnk->bfgs;
98562675beeSAlp Dener     break;
98662675beeSAlp Dener   case BNK_SCALED_GRADIENT:
98762675beeSAlp Dener     ++bnk->sgrad;
98862675beeSAlp Dener     break;
98962675beeSAlp Dener   case BNK_GRADIENT:
99062675beeSAlp Dener     ++bnk->grad;
99162675beeSAlp Dener     break;
99262675beeSAlp Dener   default:
99362675beeSAlp Dener     break;
99462675beeSAlp Dener   }
99562675beeSAlp Dener   PetscFunctionReturn(0);
99662675beeSAlp Dener }
99762675beeSAlp Dener 
99862675beeSAlp Dener /* ---------------------------------------------------------- */
99962675beeSAlp Dener 
10009b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1001eb910715SAlp Dener {
1002eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1003eb910715SAlp Dener   PetscErrorCode ierr;
10049b6ef848SAlp Dener   KSPType        ksp_type;
1005e031d6f5SAlp Dener   PetscInt       i;
1006eb910715SAlp Dener 
1007eb910715SAlp Dener   PetscFunctionBegin;
1008c4b75bccSAlp Dener   if (!tao->gradient) {
1009c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1010c4b75bccSAlp Dener   }
1011c4b75bccSAlp Dener   if (!tao->stepdirection) {
1012c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1013c4b75bccSAlp Dener   }
1014c4b75bccSAlp Dener   if (!bnk->W) {
1015c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1016c4b75bccSAlp Dener   }
1017c4b75bccSAlp Dener   if (!bnk->Xold) {
1018c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1019c4b75bccSAlp Dener   }
1020c4b75bccSAlp Dener   if (!bnk->Gold) {
1021c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1022c4b75bccSAlp Dener   }
1023c4b75bccSAlp Dener   if (!bnk->Xwork) {
1024c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1025c4b75bccSAlp Dener   }
1026c4b75bccSAlp Dener   if (!bnk->Gwork) {
1027c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1028c4b75bccSAlp Dener   }
1029c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1030c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1031c4b75bccSAlp Dener   }
1032c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1033c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1034c4b75bccSAlp Dener   }
1035c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1036c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1037c4b75bccSAlp Dener   }
1038c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1039c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1040c4b75bccSAlp Dener   }
1041e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1042c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1043c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
104489da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
104589da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
104689da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1047c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1048c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1049c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1050c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1051c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1052c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1053c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1054c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1055c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1056c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1057c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1058c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1059c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1060c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1061e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1062e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1063e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1064937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1065e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1066e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1067e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1068e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1069c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1070e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1071e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1072e031d6f5SAlp Dener     }
1073e031d6f5SAlp Dener   }
1074c0f10754SAlp Dener   bnk->X_inactive = 0;
1075c0f10754SAlp Dener   bnk->G_inactive = 0;
107662675beeSAlp Dener   bnk->inactive_work = 0;
107762675beeSAlp Dener   bnk->active_work = 0;
107862675beeSAlp Dener   bnk->inactive_idx = 0;
107962675beeSAlp Dener   bnk->active_idx = 0;
108062675beeSAlp Dener   bnk->active_lower = 0;
108162675beeSAlp Dener   bnk->active_upper = 0;
108262675beeSAlp Dener   bnk->active_fixed = 0;
1083eb910715SAlp Dener   bnk->M = 0;
108409164190SAlp Dener   bnk->H_inactive = 0;
108509164190SAlp Dener   bnk->Hpre_inactive = 0;
10869b6ef848SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
10879b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
10889b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
10899b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1090eb910715SAlp Dener   PetscFunctionReturn(0);
1091eb910715SAlp Dener }
1092eb910715SAlp Dener 
1093eb910715SAlp Dener /*------------------------------------------------------------*/
1094df278d8fSAlp Dener 
1095eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
1096eb910715SAlp Dener {
1097eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1098eb910715SAlp Dener   PetscErrorCode ierr;
1099eb910715SAlp Dener 
1100eb910715SAlp Dener   PetscFunctionBegin;
1101eb910715SAlp Dener   if (tao->setupcalled) {
1102eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1103eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1104eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
110509164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
110609164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
110709164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
110809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
110962675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
111062675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1111c4b75bccSAlp Dener   }
1112ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1113ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1114ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1115ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1116ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
1117c4b75bccSAlp Dener   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {
1118c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1119c4b75bccSAlp Dener   }
1120c4b75bccSAlp Dener   if (bnk->H_inactive != tao->hessian) {
1121c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1122c4b75bccSAlp Dener   }
1123ca964c49SAlp Dener   ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1124eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1125eb910715SAlp Dener   PetscFunctionReturn(0);
1126eb910715SAlp Dener }
1127eb910715SAlp Dener 
1128eb910715SAlp Dener /*------------------------------------------------------------*/
1129df278d8fSAlp Dener 
1130eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1131eb910715SAlp Dener {
1132eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1133eb910715SAlp Dener   PetscErrorCode ierr;
1134eb910715SAlp Dener 
1135eb910715SAlp Dener   PetscFunctionBegin;
1136eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
11378d5ead36SAlp Dener   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);
11388d5ead36SAlp Dener   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);
11392f75a4aaSAlp Dener   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);
1140748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1141748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1142748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1143748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1144748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1145748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1146748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1147748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1148748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1149748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1150748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1151748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1152748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1153748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1154748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
1155748696b1SAlp Dener   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);
1156748696b1SAlp Dener   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);
1157748696b1SAlp Dener   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);
1158748696b1SAlp Dener   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);
1159748696b1SAlp Dener   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);
1160748696b1SAlp Dener   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);
1161748696b1SAlp Dener   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);
1162748696b1SAlp Dener   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);
1163748696b1SAlp Dener   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);
1164748696b1SAlp Dener   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);
1165748696b1SAlp Dener   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);
1166748696b1SAlp Dener   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);
1167748696b1SAlp Dener   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);
1168748696b1SAlp Dener   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);
1169748696b1SAlp Dener   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);
1170748696b1SAlp Dener   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);
1171748696b1SAlp Dener   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);
1172748696b1SAlp Dener   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);
1173748696b1SAlp Dener   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);
1174748696b1SAlp Dener   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);
1175748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
1176748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1177748696b1SAlp Dener   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);
1178748696b1SAlp Dener   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);
1179748696b1SAlp Dener   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);
1180748696b1SAlp Dener   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);
1181748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1182748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1183748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1184748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1185748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
1186748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1187c0f10754SAlp Dener   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);
1188eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
118933c78596SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1190eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1191eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1192eb910715SAlp Dener   PetscFunctionReturn(0);
1193eb910715SAlp Dener }
1194eb910715SAlp Dener 
1195eb910715SAlp Dener /*------------------------------------------------------------*/
1196df278d8fSAlp Dener 
1197eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1198eb910715SAlp Dener {
1199eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1200eb910715SAlp Dener   PetscInt       nrejects;
1201eb910715SAlp Dener   PetscBool      isascii;
1202eb910715SAlp Dener   PetscErrorCode ierr;
1203eb910715SAlp Dener 
1204eb910715SAlp Dener   PetscFunctionBegin;
1205eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1206eb910715SAlp Dener   if (isascii) {
1207eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1208b9ac7092SAlp Dener     if (bnk->M) {
1209cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(bnk->M,&nrejects);CHKERRQ(ierr);
1210b9ac7092SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %D\n",nrejects);CHKERRQ(ierr);
1211eb910715SAlp Dener     }
1212e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1213eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1214b9ac7092SAlp Dener     if (bnk->M) {
1215eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1216b9ac7092SAlp Dener     }
1217eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1218eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1219eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1220eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1221eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1222eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1223eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1224eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1225eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1226eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1227eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1228eb910715SAlp Dener   }
1229eb910715SAlp Dener   PetscFunctionReturn(0);
1230eb910715SAlp Dener }
1231eb910715SAlp Dener 
1232eb910715SAlp Dener /* ---------------------------------------------------------- */
1233df278d8fSAlp Dener 
1234eb910715SAlp Dener /*MC
1235eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
123666ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1237eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1238eb910715SAlp Dener               Hk dk = -gk
12392b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12402b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1241eb910715SAlp Dener 
1242eb910715SAlp Dener     Options Database Keys:
1243*b2d8c577SAlp Dener + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
12443cb832f1SAlp Dener . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
12453cb832f1SAlp Dener . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
12463cb832f1SAlp Dener . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
1247748696b1SAlp Dener . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-tao_bnk_as_type bertsekas)
1248748696b1SAlp Dener . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-tao_bnk_as_type bertsekas)
1249748696b1SAlp Dener . -tao_bnk_sval - (developer) Hessian perturbation starting value
1250748696b1SAlp Dener . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
1251748696b1SAlp Dener . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
1252748696b1SAlp Dener . -tao_bnk_pmin - (developer) minimum Hessian perturbation
1253748696b1SAlp Dener . -tao_bnk_pmax - (developer) aximum Hessian perturbation
1254748696b1SAlp Dener . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
1255748696b1SAlp Dener . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
1256748696b1SAlp Dener . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
1257748696b1SAlp Dener . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
1258748696b1SAlp Dener . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
1259748696b1SAlp Dener . -tao_bnk_eta1 - (developer) threshold for rejecting step (-tao_bnk_update_type reduction)
1260748696b1SAlp Dener . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)
1261748696b1SAlp Dener . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)
1262748696b1SAlp Dener . -tao_bnk_eta4 - (developer) threshold for accepting good step (-tao_bnk_update_type reduction)
1263748696b1SAlp Dener . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)
1264748696b1SAlp Dener . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)
1265748696b1SAlp Dener . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)
1266748696b1SAlp Dener . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)
1267748696b1SAlp Dener . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)
1268748696b1SAlp Dener . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-tao_bnk_update_type reduction)
1269748696b1SAlp Dener . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)
1270748696b1SAlp Dener . -tao_bnk_mu2 - (developer) threshold for accepting good step (-tao_bnk_update_type interpolation)
1271748696b1SAlp Dener . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)
1272748696b1SAlp Dener . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)
1273748696b1SAlp Dener . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)
1274748696b1SAlp Dener . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)
1275748696b1SAlp Dener . -tao_bnk_theta - (developer) trust region interpolation factor (-tao_bnk_update_type interpolation)
1276748696b1SAlp Dener . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-tao_bnk_update_type step)
1277748696b1SAlp Dener . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)
1278748696b1SAlp Dener . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-tao_bnk_update_type step)
1279748696b1SAlp Dener . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-tao_bnk_update_type step)
1280748696b1SAlp Dener . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-tao_bnk_update_type step)
1281748696b1SAlp Dener . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)
1282748696b1SAlp Dener . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-tao_bnk_update_type step)
1283748696b1SAlp Dener . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)
1284748696b1SAlp Dener . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-tao_bnk_update_type step)
1285748696b1SAlp Dener . -tao_bnk_mu1_i -  (developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)
1286748696b1SAlp Dener . -tao_bnk_mu2_i -  (developer) threshold for accepting good step (-tao_bnk_init_type interpolation)
1287748696b1SAlp Dener . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-tao_bnk_init_type interpolation)
1288748696b1SAlp Dener . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-tao_bnk_init_type interpolation)
1289748696b1SAlp Dener . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-tao_bnk_init_type interpolation)
1290748696b1SAlp Dener . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-tao_bnk_init_type interpolation)
1291748696b1SAlp Dener - -tao_bnk_theta_i - (developer) trust region interpolation factor (-tao_bnk_init_type interpolation)
1292eb910715SAlp Dener 
1293eb910715SAlp Dener   Level: beginner
1294eb910715SAlp Dener M*/
1295eb910715SAlp Dener 
1296eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1297eb910715SAlp Dener {
1298eb910715SAlp Dener   TAO_BNK        *bnk;
1299eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1300eb910715SAlp Dener   PetscErrorCode ierr;
1301b9ac7092SAlp Dener   PC             pc;
1302eb910715SAlp Dener 
1303eb910715SAlp Dener   PetscFunctionBegin;
1304eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1305eb910715SAlp Dener 
1306eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1307eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1308eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1309eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1310eb910715SAlp Dener 
1311eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1312eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1313eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1314eb910715SAlp Dener 
1315eb910715SAlp Dener   tao->data = (void*)bnk;
1316eb910715SAlp Dener 
131766ed3702SAlp Dener   /*  Hessian shifting parameters */
1318eb910715SAlp Dener   bnk->sval   = 0.0;
1319eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1320eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1321eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1322eb910715SAlp Dener 
1323eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1324eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1325eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1326eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1327eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1328eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1329eb910715SAlp Dener 
1330eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1331eb910715SAlp Dener   bnk->nu1 = 0.25;
1332eb910715SAlp Dener   bnk->nu2 = 0.50;
1333eb910715SAlp Dener   bnk->nu3 = 1.00;
1334eb910715SAlp Dener   bnk->nu4 = 1.25;
1335eb910715SAlp Dener 
1336eb910715SAlp Dener   bnk->omega1 = 0.25;
1337eb910715SAlp Dener   bnk->omega2 = 0.50;
1338eb910715SAlp Dener   bnk->omega3 = 1.00;
1339eb910715SAlp Dener   bnk->omega4 = 2.00;
1340eb910715SAlp Dener   bnk->omega5 = 4.00;
1341eb910715SAlp Dener 
1342eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1343eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1344eb910715SAlp Dener   bnk->eta2 = 0.25;
1345eb910715SAlp Dener   bnk->eta3 = 0.50;
1346eb910715SAlp Dener   bnk->eta4 = 0.90;
1347eb910715SAlp Dener 
1348eb910715SAlp Dener   bnk->alpha1 = 0.25;
1349eb910715SAlp Dener   bnk->alpha2 = 0.50;
1350eb910715SAlp Dener   bnk->alpha3 = 1.00;
1351eb910715SAlp Dener   bnk->alpha4 = 2.00;
1352eb910715SAlp Dener   bnk->alpha5 = 4.00;
1353eb910715SAlp Dener 
1354eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1355eb910715SAlp Dener   bnk->mu1 = 0.10;
1356eb910715SAlp Dener   bnk->mu2 = 0.50;
1357eb910715SAlp Dener 
1358eb910715SAlp Dener   bnk->gamma1 = 0.25;
1359eb910715SAlp Dener   bnk->gamma2 = 0.50;
1360eb910715SAlp Dener   bnk->gamma3 = 2.00;
1361eb910715SAlp Dener   bnk->gamma4 = 4.00;
1362eb910715SAlp Dener 
1363eb910715SAlp Dener   bnk->theta = 0.05;
1364eb910715SAlp Dener 
1365eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1366eb910715SAlp Dener   bnk->mu1_i = 0.35;
1367eb910715SAlp Dener   bnk->mu2_i = 0.50;
1368eb910715SAlp Dener 
1369eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1370eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1371eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1372eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1373eb910715SAlp Dener 
1374eb910715SAlp Dener   bnk->theta_i = 0.25;
1375eb910715SAlp Dener 
1376eb910715SAlp Dener   /*  Remaining parameters */
1377c0f10754SAlp Dener   bnk->max_cg_its = 0;
1378eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1379eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1380770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
13810a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
13820a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
138362675beeSAlp Dener   bnk->dmin = 1.0e-6;
138462675beeSAlp Dener   bnk->dmax = 1.0e6;
1385eb910715SAlp Dener 
1386b9ac7092SAlp Dener   bnk->M               = 0;
1387b9ac7092SAlp Dener   bnk->bfgs_pre        = 0;
1388eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1389fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
13902f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1391eb910715SAlp Dener 
1392e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1393e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1394e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1395e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1396e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1397e031d6f5SAlp Dener 
1398c0f10754SAlp Dener   /* Create the line search */
1399eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1400eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1401e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1402eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1403eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1404eb910715SAlp Dener 
1405eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1406eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1407eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1408eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1409eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1410b9ac7092SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);
1411b9ac7092SAlp Dener   ierr = PCSetType(pc, PCLMVM);CHKERRQ(ierr);
1412eb910715SAlp Dener   PetscFunctionReturn(0);
1413eb910715SAlp Dener }
1414