xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 691b26d3a7ce3263bd9be9c446af0af2a46feecf)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener #include <petscksp.h>
4eb910715SAlp Dener 
570a3f44bSAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
670a3f44bSAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
770a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
870a3f44bSAlp Dener 
9e031d6f5SAlp Dener /*------------------------------------------------------------*/
10e031d6f5SAlp Dener 
11df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
12df278d8fSAlp Dener 
13c0f10754SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
14eb910715SAlp Dener {
15eb910715SAlp Dener   PetscErrorCode               ierr;
16eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
17eb910715SAlp Dener   PC                           pc;
18eb910715SAlp Dener 
1989da521bSAlp Dener   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma, resnorm;
20eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
210ad3a497SAlp Dener   PetscBool                    is_bfgs, is_jacobi, is_symmetric, sym_set;
22c4b75bccSAlp Dener   PetscInt                     n, N, nDiff;
23eb910715SAlp Dener   PetscInt                     i_max = 5;
24eb910715SAlp Dener   PetscInt                     j_max = 1;
25eb910715SAlp Dener   PetscInt                     i, j;
26eb910715SAlp Dener 
27eb910715SAlp Dener   PetscFunctionBegin;
2828017e9fSAlp Dener   /* Project the current point onto the feasible set */
2928017e9fSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
30e031d6f5SAlp Dener   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
31b9ac7092SAlp Dener   if (tao->bounded) {
3228017e9fSAlp Dener     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
33cd929ea3SAlp Dener   }
3428017e9fSAlp Dener 
3528017e9fSAlp Dener   /* Project the initial point onto the feasible region */
363b063059SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
3728017e9fSAlp Dener 
3828017e9fSAlp Dener   /* Check convergence criteria */
3928017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
4061be54a6SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
4161be54a6SAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
4261be54a6SAlp Dener   ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
43f5766c09SAlp Dener   ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
4428017e9fSAlp Dener 
45c0f10754SAlp Dener   /* Test the initial point for convergence */
4689da521bSAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
4789da521bSAlp Dener   ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
48*691b26d3SBarry Smith   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
49e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
50e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
51c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
52c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
53c0f10754SAlp Dener 
54e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
55eb910715SAlp Dener   bnk->ksp_atol = 0;
56eb910715SAlp Dener   bnk->ksp_rtol = 0;
57eb910715SAlp Dener   bnk->ksp_dtol = 0;
58eb910715SAlp Dener   bnk->ksp_ctol = 0;
59eb910715SAlp Dener   bnk->ksp_negc = 0;
60eb910715SAlp Dener   bnk->ksp_iter = 0;
61eb910715SAlp Dener   bnk->ksp_othr = 0;
62eb910715SAlp Dener 
63e031d6f5SAlp Dener   /* Reset accepted step type counters */
64e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
65e031d6f5SAlp Dener   bnk->newt = 0;
66e031d6f5SAlp Dener   bnk->bfgs = 0;
67e031d6f5SAlp Dener   bnk->sgrad = 0;
68e031d6f5SAlp Dener   bnk->grad = 0;
69e031d6f5SAlp Dener 
70fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
71fed79b8eSAlp Dener   bnk->pert = bnk->sval;
72fed79b8eSAlp Dener 
73937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
74937a31a1SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
75937a31a1SAlp Dener 
76e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
77b9ac7092SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
78b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr);
79b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr);
80b9ac7092SAlp Dener   if (is_bfgs) {
81b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
82b9ac7092SAlp Dener     ierr = PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M);CHKERRQ(ierr);
83eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
84eb910715SAlp Dener     ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
85b9ac7092SAlp Dener     ierr = MatSetSizes(bnk->M, n, n, N, N);CHKERRQ(ierr);
86cd929ea3SAlp Dener     ierr = MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
870ad3a497SAlp Dener     ierr = MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric);CHKERRQ(ierr);
880ad3a497SAlp Dener     if (!sym_set || !is_symmetric) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
89b9ac7092SAlp Dener   } else if (is_jacobi) {
90b9ac7092SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
91e031d6f5SAlp Dener   }
92e031d6f5SAlp Dener 
93e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
9462675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
9562675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
96eb910715SAlp Dener 
97eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
98eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
99c0f10754SAlp Dener   *needH = PETSC_TRUE;
100eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
10162675beeSAlp Dener     switch(initType) {
102eb910715SAlp Dener     case BNK_INIT_CONSTANT:
103eb910715SAlp Dener       /* Use the initial radius specified */
104c0f10754SAlp Dener       tao->trust = tao->trust0;
105eb910715SAlp Dener       break;
106eb910715SAlp Dener 
107eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
108c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
109eb910715SAlp Dener       max_radius = 0.0;
11008752603SAlp Dener       tao->trust = tao->trust0;
111eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1120a4511e9SAlp Dener         f_min = bnk->f;
113eb910715SAlp Dener         sigma = 0.0;
114eb910715SAlp Dener 
115c0f10754SAlp Dener         if (*needH) {
11662602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
117e0ed867bSAlp Dener           ierr = bnk->computehessian(tao);CHKERRQ(ierr);
11808752603SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
11989da521bSAlp Dener           ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
12089da521bSAlp Dener           if (bnk->active_idx) {
1212ab2a32cSAlp Dener             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
12228017e9fSAlp Dener           } else {
12308752603SAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
12428017e9fSAlp Dener           }
125c0f10754SAlp Dener           *needH = PETSC_FALSE;
126eb910715SAlp Dener         }
127eb910715SAlp Dener 
128eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
12962602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
13062602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
13162602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
1323b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
13389da521bSAlp Dener           /* Compute the step we actually accepted */
134eb910715SAlp Dener           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
13562602cfbSAlp Dener           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
13662602cfbSAlp Dener           /* Compute the objective at the trial */
13762602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
138*691b26d3SBarry Smith           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
13962602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
140eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
141eb910715SAlp Dener             tau = bnk->gamma1_i;
142eb910715SAlp Dener           } else {
1430a4511e9SAlp Dener             if (ftrial < f_min) {
1440a4511e9SAlp Dener               f_min = ftrial;
145eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
146eb910715SAlp Dener             }
14708752603SAlp Dener 
148770b7498SAlp Dener             /* Compute the predicted and actual reduction */
14989da521bSAlp Dener             if (bnk->active_idx) {
15008752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
15108752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1522ab2a32cSAlp Dener             } else {
15308752603SAlp Dener               bnk->X_inactive = bnk->W;
15408752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1552ab2a32cSAlp Dener             }
15608752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
15708752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
15889da521bSAlp Dener             if (bnk->active_idx) {
15908752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
16008752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1612ab2a32cSAlp Dener             }
162eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
163eb910715SAlp Dener             actred = bnk->f - ftrial;
1643105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
165eb910715SAlp Dener               kappa = 1.0;
1663105154fSTodd Munson             } else {
167eb910715SAlp Dener               kappa = actred / prered;
168eb910715SAlp Dener             }
169eb910715SAlp Dener 
170eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
171eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
172eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
173eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
174eb910715SAlp Dener 
17518cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
176eb910715SAlp Dener               /*  Great agreement */
177eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
178eb910715SAlp Dener 
179eb910715SAlp Dener               if (tau_max < 1.0) {
180eb910715SAlp Dener                 tau = bnk->gamma3_i;
1813105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
182eb910715SAlp Dener                 tau = bnk->gamma4_i;
1833105154fSTodd Munson               } else {
184eb910715SAlp Dener                 tau = tau_max;
185eb910715SAlp Dener               }
18618cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
187eb910715SAlp Dener               /*  Good agreement */
188eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
189eb910715SAlp Dener 
190eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
191eb910715SAlp Dener                 tau = bnk->gamma2_i;
192eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
193eb910715SAlp Dener                 tau = bnk->gamma3_i;
194eb910715SAlp Dener               } else {
195eb910715SAlp Dener                 tau = tau_max;
196eb910715SAlp Dener               }
1978f8a4e06SAlp Dener             } else {
198eb910715SAlp Dener               /*  Not good agreement */
199eb910715SAlp Dener               if (tau_min > 1.0) {
200eb910715SAlp Dener                 tau = bnk->gamma2_i;
201eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
202eb910715SAlp Dener                 tau = bnk->gamma1_i;
203eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
204eb910715SAlp Dener                 tau = bnk->gamma1_i;
2053105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
206eb910715SAlp Dener                 tau = tau_1;
2073105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
208eb910715SAlp Dener                 tau = tau_2;
209eb910715SAlp Dener               } else {
210eb910715SAlp Dener                 tau = tau_max;
211eb910715SAlp Dener               }
212eb910715SAlp Dener             }
213eb910715SAlp Dener           }
214eb910715SAlp Dener           tao->trust = tau * tao->trust;
215eb910715SAlp Dener         }
216eb910715SAlp Dener 
2170a4511e9SAlp Dener         if (f_min < bnk->f) {
218937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2190a4511e9SAlp Dener           bnk->f = f_min;
220937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
221eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
2223b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
223937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
224937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
22509164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
22661be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
22761be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
22861be54a6SAlp Dener           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
229937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
230f5766c09SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
231c0f10754SAlp Dener           *needH = PETSC_TRUE;
232937a31a1SAlp Dener           /* Test the new step for convergence */
23389da521bSAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
23489da521bSAlp Dener           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
235*691b26d3SBarry Smith           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
236e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
237e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
238eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
239eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
240937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
241937a31a1SAlp Dener           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
242eb910715SAlp Dener         }
243eb910715SAlp Dener       }
244eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
245e031d6f5SAlp Dener 
246e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
247e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
248e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
249eb910715SAlp Dener       break;
250eb910715SAlp Dener 
251eb910715SAlp Dener     default:
252eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
253eb910715SAlp Dener       tao->trust = 0.0;
254eb910715SAlp Dener       break;
255eb910715SAlp Dener     }
256eb910715SAlp Dener   }
257eb910715SAlp Dener   PetscFunctionReturn(0);
258eb910715SAlp Dener }
259eb910715SAlp Dener 
260df278d8fSAlp Dener /*------------------------------------------------------------*/
261df278d8fSAlp Dener 
262e0ed867bSAlp Dener /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */
26362675beeSAlp Dener 
26462675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
26562675beeSAlp Dener {
26662675beeSAlp Dener   PetscErrorCode               ierr;
26762675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
26862675beeSAlp Dener 
26962675beeSAlp Dener   PetscFunctionBegin;
27062675beeSAlp Dener   /* Compute the Hessian */
27162675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
27262675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
273b9ac7092SAlp Dener   if (bnk->M) {
27462675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
27562675beeSAlp Dener   }
276e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
277f5766c09SAlp Dener   if (bnk->Hpre_inactive) {
278f5766c09SAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
279f5766c09SAlp Dener   }
280f5766c09SAlp Dener   if (bnk->H_inactive) {
281e0ed867bSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
282f5766c09SAlp Dener   }
283f5766c09SAlp Dener   if (bnk->active_idx) {
284e0ed867bSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
285e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
286f5766c09SAlp Dener       ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
287e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
288e0ed867bSAlp Dener     } else {
289e0ed867bSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
290e0ed867bSAlp Dener     }
291e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
292e0ed867bSAlp Dener       ierr = PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx);CHKERRQ(ierr);
293e0ed867bSAlp Dener     }
294e0ed867bSAlp Dener   } else {
295e0ed867bSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
296e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
297f5766c09SAlp Dener       ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
298e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
299e0ed867bSAlp Dener     } else {
300e0ed867bSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);CHKERRQ(ierr);
301e0ed867bSAlp Dener     }
302e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
303e0ed867bSAlp Dener       ierr = PCLMVMClearIS(bnk->bfgs_pre);CHKERRQ(ierr);
304e0ed867bSAlp Dener     }
305e0ed867bSAlp Dener   }
30662675beeSAlp Dener   PetscFunctionReturn(0);
30762675beeSAlp Dener }
30862675beeSAlp Dener 
30962675beeSAlp Dener /*------------------------------------------------------------*/
31062675beeSAlp Dener 
3112f75a4aaSAlp Dener /* Routine for estimating the active set */
3122f75a4aaSAlp Dener 
31308752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3142f75a4aaSAlp Dener {
3152f75a4aaSAlp Dener   PetscErrorCode               ierr;
3162f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
317c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3182f75a4aaSAlp Dener 
3192f75a4aaSAlp Dener   PetscFunctionBegin;
32008752603SAlp Dener   switch (asType) {
3212f75a4aaSAlp Dener   case BNK_AS_NONE:
3222f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3232f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3242f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3252f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3262f75a4aaSAlp Dener     break;
3272f75a4aaSAlp Dener 
3282f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3292f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
330b9ac7092SAlp Dener     if (bnk->M) {
3312f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3329515a401SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3332f75a4aaSAlp Dener     } else {
334fc5ca067SStefano Zampini       hessComputed = diagExists = PETSC_FALSE;
335f5766c09SAlp Dener       if (tao->hessian) {
33661be54a6SAlp Dener         ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
337f5766c09SAlp Dener       }
338fc5ca067SStefano Zampini       if (hessComputed) {
339fc5ca067SStefano Zampini         ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
340fc5ca067SStefano Zampini       }
341fc5ca067SStefano Zampini       if (diagExists) {
3429b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3432f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3449b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3459b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3462f75a4aaSAlp Dener         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3472f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
34861be54a6SAlp Dener       } else {
349c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
35061be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
35161be54a6SAlp Dener       }
3522f75a4aaSAlp Dener     }
3532f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
35489da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
35589da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
356c4b75bccSAlp Dener     break;
3572f75a4aaSAlp Dener 
3582f75a4aaSAlp Dener   default:
3592f75a4aaSAlp Dener     break;
3602f75a4aaSAlp Dener   }
3612f75a4aaSAlp Dener   PetscFunctionReturn(0);
3622f75a4aaSAlp Dener }
3632f75a4aaSAlp Dener 
3642f75a4aaSAlp Dener /*------------------------------------------------------------*/
3652f75a4aaSAlp Dener 
3662f75a4aaSAlp Dener /* Routine for bounding the step direction */
3672f75a4aaSAlp Dener 
368a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
3692f75a4aaSAlp Dener {
3702f75a4aaSAlp Dener   PetscErrorCode               ierr;
3712f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3722f75a4aaSAlp Dener 
3732f75a4aaSAlp Dener   PetscFunctionBegin;
374a1318120SAlp Dener   switch (asType) {
3752f75a4aaSAlp Dener   case BNK_AS_NONE:
376c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
3772f75a4aaSAlp Dener     break;
3782f75a4aaSAlp Dener 
3792f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
380c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
3812f75a4aaSAlp Dener     break;
3822f75a4aaSAlp Dener 
3832f75a4aaSAlp Dener   default:
3842f75a4aaSAlp Dener     break;
3852f75a4aaSAlp Dener   }
3862f75a4aaSAlp Dener   PetscFunctionReturn(0);
3872f75a4aaSAlp Dener }
3882f75a4aaSAlp Dener 
389e031d6f5SAlp Dener /*------------------------------------------------------------*/
390e031d6f5SAlp Dener 
391e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
392e031d6f5SAlp Dener    accelerate Newton convergence.
393e031d6f5SAlp Dener 
394e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
395e031d6f5SAlp Dener    for more gradient evaluations.
396e031d6f5SAlp Dener */
397e031d6f5SAlp Dener 
398c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
399c0f10754SAlp Dener {
400c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
401c0f10754SAlp Dener   PetscErrorCode               ierr;
402c0f10754SAlp Dener 
403c0f10754SAlp Dener   PetscFunctionBegin;
404c0f10754SAlp Dener   *terminate = PETSC_FALSE;
405c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
406c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
407c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
408c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
409c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
410c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
411c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
412c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
413c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
414c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
415e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
416c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
417c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
418c0f10754SAlp 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) {
419c0f10754SAlp Dener       *terminate = PETSC_TRUE;
42061be54a6SAlp Dener     } else {
42133c78596SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
422c0f10754SAlp Dener     }
423c0f10754SAlp Dener   }
424c0f10754SAlp Dener   PetscFunctionReturn(0);
425c0f10754SAlp Dener }
426c0f10754SAlp Dener 
4272f75a4aaSAlp Dener /*------------------------------------------------------------*/
4282f75a4aaSAlp Dener 
429c0f10754SAlp Dener /* Routine for computing the Newton step. */
430df278d8fSAlp Dener 
4316b591159SAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
432eb910715SAlp Dener {
433eb910715SAlp Dener   PetscErrorCode               ierr;
434eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
435eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
436eb910715SAlp Dener   PetscInt                     kspits;
437bddd1ffdSAlp Dener   PetscBool                    is_lmvm;
438eb910715SAlp Dener 
439eb910715SAlp Dener   PetscFunctionBegin;
44089da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
44189da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
44289da521bSAlp Dener   if (!bnk->inactive_idx) {
44389da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
444a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
44589da521bSAlp Dener     PetscFunctionReturn(0);
44689da521bSAlp Dener   }
44789da521bSAlp Dener 
44862675beeSAlp Dener   /* Shift the reduced Hessian matrix */
44962675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
450f7bf01afSAlp Dener     ierr = PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm);CHKERRQ(ierr);
451f7bf01afSAlp Dener     if (is_lmvm) {
452f7bf01afSAlp Dener       ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
453f7bf01afSAlp Dener     } else {
45462675beeSAlp Dener       ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
45562675beeSAlp Dener       if (bnk->H_inactive != bnk->Hpre_inactive) {
45662675beeSAlp Dener         ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
45762675beeSAlp Dener       }
45862675beeSAlp Dener     }
459f7bf01afSAlp Dener   }
46062675beeSAlp Dener 
461eb910715SAlp Dener   /* Solve the Newton system of equations */
462937a31a1SAlp Dener   tao->ksp_its = 0;
4632f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
4645e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
46509164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
4665e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
46789da521bSAlp Dener   if (bnk->active_idx) {
4685e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
4695e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
4705e9b73cbSAlp Dener   } else {
4715e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4725e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
47328017e9fSAlp Dener   }
474eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
475fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4765e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
477eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
478eb910715SAlp Dener     tao->ksp_its+=kspits;
479eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
480080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
481eb910715SAlp Dener 
482eb910715SAlp Dener     if (0.0 == tao->trust) {
483eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
484080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
485080d2917SAlp Dener         tao->trust = bnk->dnorm;
486eb910715SAlp Dener 
487eb910715SAlp Dener         /* Modify the radius if it is too large or small */
488eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
489eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
490eb910715SAlp Dener       } else {
491eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
492eb910715SAlp Dener            the trust-region subproblem to get a direction */
493eb910715SAlp Dener         tao->trust = tao->trust0;
494eb910715SAlp Dener 
495eb910715SAlp Dener         /* Modify the radius if it is too large or small */
496eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
497eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
498eb910715SAlp Dener 
499fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5005e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
501eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
502eb910715SAlp Dener         tao->ksp_its+=kspits;
503eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
504080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
505eb910715SAlp Dener 
506*691b26d3SBarry Smith         if (bnk->dnorm == 0.0) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero");
507eb910715SAlp Dener       }
508eb910715SAlp Dener     }
509eb910715SAlp Dener   } else {
5105e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
511eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
512eb910715SAlp Dener     tao->ksp_its += kspits;
513eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
514eb910715SAlp Dener   }
5155e9b73cbSAlp Dener   /* Restore sub vectors back */
51689da521bSAlp Dener   if (bnk->active_idx) {
5175e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5185e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5195e9b73cbSAlp Dener   }
520770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
521fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
522a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
523770b7498SAlp Dener 
524770b7498SAlp Dener   /* Record convergence reasons */
525e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
526e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
527770b7498SAlp Dener     ++bnk->ksp_atol;
528e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
529770b7498SAlp Dener     ++bnk->ksp_rtol;
530e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
531770b7498SAlp Dener     ++bnk->ksp_ctol;
532e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
533770b7498SAlp Dener     ++bnk->ksp_negc;
534e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
535770b7498SAlp Dener     ++bnk->ksp_dtol;
536e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
537770b7498SAlp Dener     ++bnk->ksp_iter;
538770b7498SAlp Dener   } else {
539770b7498SAlp Dener     ++bnk->ksp_othr;
540770b7498SAlp Dener   }
541fed79b8eSAlp Dener 
542fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
543b9ac7092SAlp Dener   if (bnk->M) {
544cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
545b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
546fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
547cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
54809164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
549eb910715SAlp Dener     }
550fed79b8eSAlp Dener   }
5516b591159SAlp Dener   *step_type = BNK_NEWTON;
552e465cd6fSAlp Dener   PetscFunctionReturn(0);
553e465cd6fSAlp Dener }
554eb910715SAlp Dener 
55562675beeSAlp Dener /*------------------------------------------------------------*/
55662675beeSAlp Dener 
5575e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5585e9b73cbSAlp Dener 
5595e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
5605e9b73cbSAlp Dener {
5615e9b73cbSAlp Dener   PetscErrorCode               ierr;
5625e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
5635e9b73cbSAlp Dener 
5645e9b73cbSAlp Dener   PetscFunctionBegin;
5655e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
56689da521bSAlp Dener   if (bnk->active_idx){
5675e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5685e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5695e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5705e9b73cbSAlp Dener   } else {
5715e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
5725e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5735e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
5745e9b73cbSAlp Dener   }
5755e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5765e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
5775e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
57833c78596SAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);CHKERRQ(ierr);
5795e9b73cbSAlp Dener   /* Restore the sub vectors */
58089da521bSAlp Dener   if (bnk->active_idx){
5815e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5825e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5835e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5845e9b73cbSAlp Dener   }
5855e9b73cbSAlp Dener   PetscFunctionReturn(0);
5865e9b73cbSAlp Dener }
5875e9b73cbSAlp Dener 
5885e9b73cbSAlp Dener /*------------------------------------------------------------*/
5895e9b73cbSAlp Dener 
59062675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
59162675beeSAlp Dener 
59262675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
59362675beeSAlp Dener    in the event that the Newton step fails the test.
59462675beeSAlp Dener */
59562675beeSAlp Dener 
596e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
597e465cd6fSAlp Dener {
598e465cd6fSAlp Dener   PetscErrorCode               ierr;
599e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
600e465cd6fSAlp Dener 
601b2d8c577SAlp Dener   PetscReal                    gdx, e_min;
602e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
603e465cd6fSAlp Dener 
604e465cd6fSAlp Dener   PetscFunctionBegin;
6056b591159SAlp Dener   switch (*stepType) {
6066b591159SAlp Dener   case BNK_NEWTON:
607080d2917SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
608eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
609eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
610eb910715SAlp Dener         Update the perturbation for next time */
611eb910715SAlp Dener       if (bnk->pert <= 0.0) {
612eb910715SAlp Dener         /* Initialize the perturbation */
613eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
614eb910715SAlp Dener         if (bnk->is_gltr) {
61505de396fSBarry Smith           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
616eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
617eb910715SAlp Dener         }
618eb910715SAlp Dener       } else {
619eb910715SAlp Dener         /* Increase the perturbation */
620eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
621eb910715SAlp Dener       }
622eb910715SAlp Dener 
6230ad3a497SAlp Dener       if (!bnk->M) {
624eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
625eb910715SAlp Dener           Must use gradient direction in this case */
626080d2917SAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
627eb910715SAlp Dener         *stepType = BNK_GRADIENT;
628eb910715SAlp Dener       } else {
629eb910715SAlp Dener         /* Attempt to use the BFGS direction */
6309515a401SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
631eb910715SAlp Dener 
6328d5ead36SAlp Dener         /* Check for success (descent direction)
6338d5ead36SAlp Dener           NOTE: Negative gdx here means not a descent direction because
6348d5ead36SAlp Dener           the fall-back step is missing a negative sign. */
635080d2917SAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
6363105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
637eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
638eb910715SAlp Dener             We can assert bfgsUpdates > 1 in this case because
639eb910715SAlp Dener             the first solve produces the scaled gradient direction,
640eb910715SAlp Dener             which is guaranteed to be descent */
641eb910715SAlp Dener 
642eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
643cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
64409164190SAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
6459515a401SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
646eb910715SAlp Dener 
647eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
648eb910715SAlp Dener         } else {
649cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
650eb910715SAlp Dener           if (1 == bfgsUpdates) {
651eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
652eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
653eb910715SAlp Dener           } else {
654eb910715SAlp Dener             *stepType = BNK_BFGS;
655eb910715SAlp Dener           }
656eb910715SAlp Dener         }
657eb910715SAlp Dener       }
6588d5ead36SAlp Dener       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6598d5ead36SAlp Dener       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
660a1318120SAlp Dener       ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
661eb910715SAlp Dener     } else {
662eb910715SAlp Dener       /* Computed Newton step is descent */
663eb910715SAlp Dener       switch (ksp_reason) {
664eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
665eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
666eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
667eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
668eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
669eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
670eb910715SAlp Dener         if (bnk->pert <= 0.0) {
671eb910715SAlp Dener           /* Initialize the perturbation */
672eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
673eb910715SAlp Dener           if (bnk->is_gltr) {
67405de396fSBarry Smith             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
675eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
676eb910715SAlp Dener           }
677eb910715SAlp Dener         } else {
678eb910715SAlp Dener           /* Increase the perturbation */
679eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
680eb910715SAlp Dener         }
681eb910715SAlp Dener         break;
682eb910715SAlp Dener 
683eb910715SAlp Dener       default:
684eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
685eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
686eb910715SAlp Dener         if (bnk->pert < bnk->pmin) {
687eb910715SAlp Dener           bnk->pert = 0.0;
688eb910715SAlp Dener         }
689eb910715SAlp Dener         break;
690eb910715SAlp Dener       }
691fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
692eb910715SAlp Dener     }
6936b591159SAlp Dener     break;
6946b591159SAlp Dener 
6956b591159SAlp Dener   case BNK_BFGS:
6966b591159SAlp Dener     /* Check for success (descent direction) */
6976b591159SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
6986b591159SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
6996b591159SAlp Dener       /* Step is not descent or solve was not successful
7006b591159SAlp Dener          Use steepest descent direction (scaled) */
7016b591159SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
7026b591159SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
7039515a401SAlp Dener       ierr = MatSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
7046b591159SAlp Dener       ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
7056b591159SAlp Dener       ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
7066b591159SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
7076b591159SAlp Dener     } else {
7086b591159SAlp Dener       *stepType = BNK_BFGS;
7096b591159SAlp Dener     }
7106b591159SAlp Dener     break;
7116b591159SAlp Dener 
7126b591159SAlp Dener   case BNK_SCALED_GRADIENT:
7136b591159SAlp Dener     break;
7146b591159SAlp Dener 
7156b591159SAlp Dener   default:
7166b591159SAlp Dener     break;
7176b591159SAlp Dener   }
7186b591159SAlp Dener 
719eb910715SAlp Dener   PetscFunctionReturn(0);
720eb910715SAlp Dener }
721eb910715SAlp Dener 
722df278d8fSAlp Dener /*------------------------------------------------------------*/
723df278d8fSAlp Dener 
724df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
725df278d8fSAlp Dener 
726df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
727df278d8fSAlp Dener   Newton step does not produce a valid step length.
728df278d8fSAlp Dener */
729df278d8fSAlp Dener 
730937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
731c14b763aSAlp Dener {
732c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
733c14b763aSAlp Dener   PetscErrorCode ierr;
734c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
735c14b763aSAlp Dener 
736b2d8c577SAlp Dener   PetscReal      e_min, gdx;
737c14b763aSAlp Dener   PetscInt       bfgsUpdates;
738c14b763aSAlp Dener 
739c14b763aSAlp Dener   PetscFunctionBegin;
740c14b763aSAlp Dener   /* Perform the linesearch */
741c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
742c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
743c14b763aSAlp Dener 
744b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
745c14b763aSAlp Dener     /* Linesearch failed, revert solution */
746c14b763aSAlp Dener     bnk->f = bnk->fold;
747c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
748c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
749c14b763aSAlp Dener 
750937a31a1SAlp Dener     switch(*stepType) {
751c14b763aSAlp Dener     case BNK_NEWTON:
7528d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
753c14b763aSAlp Dener          Update the perturbation for next time */
754c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
755c14b763aSAlp Dener         /* Initialize the perturbation */
756c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
757c14b763aSAlp Dener         if (bnk->is_gltr) {
75805de396fSBarry Smith           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
759c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
760c14b763aSAlp Dener         }
761c14b763aSAlp Dener       } else {
762c14b763aSAlp Dener         /* Increase the perturbation */
763c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
764c14b763aSAlp Dener       }
765c14b763aSAlp Dener 
7660ad3a497SAlp Dener       if (!bnk->M) {
767c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
768c14b763aSAlp Dener            Must use gradient direction in this case */
769937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
770937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
771c14b763aSAlp Dener       } else {
772c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
7739515a401SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
7748d5ead36SAlp Dener         /* Check for success (descent direction)
7758d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
776c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7773105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
778c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
779c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
780c14b763aSAlp Dener              Use steepest descent direction (scaled) */
781cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
782c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
7839515a401SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
784c14b763aSAlp Dener 
785c14b763aSAlp Dener           bfgsUpdates = 1;
786937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
787c14b763aSAlp Dener         } else {
788cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
789c14b763aSAlp Dener           if (1 == bfgsUpdates) {
790c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
791937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
792c14b763aSAlp Dener           } else {
793937a31a1SAlp Dener             *stepType = BNK_BFGS;
794c14b763aSAlp Dener           }
795c14b763aSAlp Dener         }
796c14b763aSAlp Dener       }
797c14b763aSAlp Dener       break;
798c14b763aSAlp Dener 
799c14b763aSAlp Dener     case BNK_BFGS:
800c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
801c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
802c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
803cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
804c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
8059515a401SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
806c14b763aSAlp Dener 
807c14b763aSAlp Dener       bfgsUpdates = 1;
808937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
809c14b763aSAlp Dener       break;
810c14b763aSAlp Dener     }
8118d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8128d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
813a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
814c14b763aSAlp Dener 
8158d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
816c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
817c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
818c14b763aSAlp Dener   }
819c14b763aSAlp Dener   *reason = ls_reason;
820c14b763aSAlp Dener   PetscFunctionReturn(0);
821c14b763aSAlp Dener }
822c14b763aSAlp Dener 
823df278d8fSAlp Dener /*------------------------------------------------------------*/
824df278d8fSAlp Dener 
825df278d8fSAlp Dener /* Routine for updating the trust radius.
826df278d8fSAlp Dener 
827df278d8fSAlp Dener   Function features three different update methods:
828df278d8fSAlp Dener   1) Line-search step length based
829df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
830df278d8fSAlp Dener   3) Interpolation
831df278d8fSAlp Dener */
832df278d8fSAlp Dener 
83328017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
834080d2917SAlp Dener {
835080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
836080d2917SAlp Dener   PetscErrorCode ierr;
837080d2917SAlp Dener 
838b1c2d0e3SAlp Dener   PetscReal      step, kappa;
839080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
840080d2917SAlp Dener 
841080d2917SAlp Dener   PetscFunctionBegin;
842080d2917SAlp Dener   /* Update trust region radius */
843080d2917SAlp Dener   *accept = PETSC_FALSE;
84428017e9fSAlp Dener   switch(updateType) {
845080d2917SAlp Dener   case BNK_UPDATE_STEP:
846c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
847080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
848080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
849080d2917SAlp Dener       if (step < bnk->nu1) {
850080d2917SAlp Dener         /* Very bad step taken; reduce radius */
851080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
852080d2917SAlp Dener       } else if (step < bnk->nu2) {
853080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
854080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
855080d2917SAlp Dener       } else if (step < bnk->nu3) {
856080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
857080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
858080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
859080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
860080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
861080d2917SAlp Dener         }
862080d2917SAlp Dener       } else if (step < bnk->nu4) {
863080d2917SAlp Dener         /*  Full step taken; increase the radius */
864080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
865080d2917SAlp Dener       } else {
866080d2917SAlp Dener         /*  More than full step taken; increase the radius */
867080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
868080d2917SAlp Dener       }
869080d2917SAlp Dener     } else {
870080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
871080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
872080d2917SAlp Dener     }
873080d2917SAlp Dener     break;
874080d2917SAlp Dener 
875080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
876080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
877e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
878fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
879fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
880fed79b8eSAlp Dener            be rejected! */
881080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8823105154fSTodd Munson       } else {
883b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
884080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
885080d2917SAlp Dener         } else {
8863105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
887080d2917SAlp Dener             kappa = 1.0;
8883105154fSTodd Munson           } else {
889080d2917SAlp Dener             kappa = actred / prered;
890080d2917SAlp Dener           }
891fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
892080d2917SAlp Dener           if (kappa < bnk->eta1) {
893fed79b8eSAlp Dener             /* Reject the step */
894080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8953105154fSTodd Munson           } else {
896fed79b8eSAlp Dener             /* Accept the step */
897c133c014SAlp Dener             *accept = PETSC_TRUE;
898c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8998d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
900080d2917SAlp Dener               if (kappa < bnk->eta2) {
901080d2917SAlp Dener                 /* Marginal bad step */
902c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
9033105154fSTodd Munson               } else if (kappa < bnk->eta3) {
904fed79b8eSAlp Dener                 /* Reasonable step */
905fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
9063105154fSTodd Munson               } else if (kappa < bnk->eta4) {
907080d2917SAlp Dener                 /* Good step */
908c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
9093105154fSTodd Munson               } else {
910080d2917SAlp Dener                 /* Very good step */
911c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
912080d2917SAlp Dener               }
913c133c014SAlp Dener             }
914080d2917SAlp Dener           }
915080d2917SAlp Dener         }
916080d2917SAlp Dener       }
917080d2917SAlp Dener     } else {
918080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
919080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
920080d2917SAlp Dener     }
921080d2917SAlp Dener     break;
922080d2917SAlp Dener 
923080d2917SAlp Dener   default:
924080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
925b1c2d0e3SAlp Dener       if (prered < 0.0) {
926080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
927080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
928080d2917SAlp Dener         /*  be rejected! */
929080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
930080d2917SAlp Dener       } else {
931b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
932080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
933080d2917SAlp Dener         } else {
934080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
935080d2917SAlp Dener             kappa = 1.0;
936080d2917SAlp Dener           } else {
937080d2917SAlp Dener             kappa = actred / prered;
938080d2917SAlp Dener           }
939080d2917SAlp Dener 
940080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
941080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
942080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
943080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
944080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
945080d2917SAlp Dener 
946080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
947080d2917SAlp Dener             /*  Great agreement */
948080d2917SAlp Dener             *accept = PETSC_TRUE;
949080d2917SAlp Dener             if (tau_max < 1.0) {
950080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
951080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
952080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
953080d2917SAlp Dener             } else {
954080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
955080d2917SAlp Dener             }
956080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
957080d2917SAlp Dener             /*  Good agreement */
958080d2917SAlp Dener             *accept = PETSC_TRUE;
959080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
960080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
961080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
962080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
963080d2917SAlp Dener             } else if (tau_max < 1.0) {
964080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
965080d2917SAlp Dener             } else {
966080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
967080d2917SAlp Dener             }
968080d2917SAlp Dener           } else {
969080d2917SAlp Dener             /*  Not good agreement */
970080d2917SAlp Dener             if (tau_min > 1.0) {
971080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
972080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
973080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
974080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
975080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
976080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
977080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
978080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
979080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
980080d2917SAlp Dener             } else {
981080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
982080d2917SAlp Dener             }
983080d2917SAlp Dener           }
984080d2917SAlp Dener         }
985080d2917SAlp Dener       }
986080d2917SAlp Dener     } else {
987080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
988080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
989080d2917SAlp Dener     }
99028017e9fSAlp Dener     break;
991080d2917SAlp Dener   }
992c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
993c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
994fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
995080d2917SAlp Dener   PetscFunctionReturn(0);
996080d2917SAlp Dener }
997080d2917SAlp Dener 
998eb910715SAlp Dener /* ---------------------------------------------------------- */
999df278d8fSAlp Dener 
100062675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
100162675beeSAlp Dener {
100262675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
100362675beeSAlp Dener 
100462675beeSAlp Dener   PetscFunctionBegin;
100562675beeSAlp Dener   switch (stepType) {
100662675beeSAlp Dener   case BNK_NEWTON:
100762675beeSAlp Dener     ++bnk->newt;
100862675beeSAlp Dener     break;
100962675beeSAlp Dener   case BNK_BFGS:
101062675beeSAlp Dener     ++bnk->bfgs;
101162675beeSAlp Dener     break;
101262675beeSAlp Dener   case BNK_SCALED_GRADIENT:
101362675beeSAlp Dener     ++bnk->sgrad;
101462675beeSAlp Dener     break;
101562675beeSAlp Dener   case BNK_GRADIENT:
101662675beeSAlp Dener     ++bnk->grad;
101762675beeSAlp Dener     break;
101862675beeSAlp Dener   default:
101962675beeSAlp Dener     break;
102062675beeSAlp Dener   }
102162675beeSAlp Dener   PetscFunctionReturn(0);
102262675beeSAlp Dener }
102362675beeSAlp Dener 
102462675beeSAlp Dener /* ---------------------------------------------------------- */
102562675beeSAlp Dener 
10269b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1027eb910715SAlp Dener {
1028eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1029eb910715SAlp Dener   PetscErrorCode ierr;
1030e031d6f5SAlp Dener   PetscInt       i;
1031eb910715SAlp Dener 
1032eb910715SAlp Dener   PetscFunctionBegin;
1033c4b75bccSAlp Dener   if (!tao->gradient) {
1034c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1035c4b75bccSAlp Dener   }
1036c4b75bccSAlp Dener   if (!tao->stepdirection) {
1037c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1038c4b75bccSAlp Dener   }
1039c4b75bccSAlp Dener   if (!bnk->W) {
1040c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1041c4b75bccSAlp Dener   }
1042c4b75bccSAlp Dener   if (!bnk->Xold) {
1043c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1044c4b75bccSAlp Dener   }
1045c4b75bccSAlp Dener   if (!bnk->Gold) {
1046c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1047c4b75bccSAlp Dener   }
1048c4b75bccSAlp Dener   if (!bnk->Xwork) {
1049c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1050c4b75bccSAlp Dener   }
1051c4b75bccSAlp Dener   if (!bnk->Gwork) {
1052c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1053c4b75bccSAlp Dener   }
1054c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1055c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1056c4b75bccSAlp Dener   }
1057c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1058c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1059c4b75bccSAlp Dener   }
1060c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1061c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1062c4b75bccSAlp Dener   }
1063c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1064c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1065c4b75bccSAlp Dener   }
1066e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1067c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1068c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
106989da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
107089da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
107189da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1072c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1073c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1074c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1075c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1076c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1077c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1078c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1079c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1080c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1081c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1082c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1083c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1084c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1085c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1086e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1087e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1088e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1089937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1090e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1091e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1092e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1093e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1094c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1095e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1096e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1097e031d6f5SAlp Dener     }
1098e031d6f5SAlp Dener   }
1099c0f10754SAlp Dener   bnk->X_inactive = 0;
1100c0f10754SAlp Dener   bnk->G_inactive = 0;
110162675beeSAlp Dener   bnk->inactive_work = 0;
110262675beeSAlp Dener   bnk->active_work = 0;
110362675beeSAlp Dener   bnk->inactive_idx = 0;
110462675beeSAlp Dener   bnk->active_idx = 0;
110562675beeSAlp Dener   bnk->active_lower = 0;
110662675beeSAlp Dener   bnk->active_upper = 0;
110762675beeSAlp Dener   bnk->active_fixed = 0;
1108eb910715SAlp Dener   bnk->M = 0;
110909164190SAlp Dener   bnk->H_inactive = 0;
111009164190SAlp Dener   bnk->Hpre_inactive = 0;
1111eb910715SAlp Dener   PetscFunctionReturn(0);
1112eb910715SAlp Dener }
1113eb910715SAlp Dener 
1114eb910715SAlp Dener /*------------------------------------------------------------*/
1115df278d8fSAlp Dener 
1116e0ed867bSAlp Dener PetscErrorCode TaoDestroy_BNK(Tao tao)
1117eb910715SAlp Dener {
1118eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1119eb910715SAlp Dener   PetscErrorCode ierr;
1120eb910715SAlp Dener 
1121eb910715SAlp Dener   PetscFunctionBegin;
1122eb910715SAlp Dener   if (tao->setupcalled) {
1123eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1124eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1125eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
112609164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
112709164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
112809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
112909164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
113062675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
113162675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1132c4b75bccSAlp Dener   }
1133ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1134ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1135ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1136ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1137ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
1138c4b75bccSAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1139c4b75bccSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1140ca964c49SAlp Dener   ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1141eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1142eb910715SAlp Dener   PetscFunctionReturn(0);
1143eb910715SAlp Dener }
1144eb910715SAlp Dener 
1145eb910715SAlp Dener /*------------------------------------------------------------*/
1146df278d8fSAlp Dener 
1147e0ed867bSAlp Dener PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1148eb910715SAlp Dener {
1149eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1150eb910715SAlp Dener   PetscErrorCode ierr;
1151e0ed867bSAlp Dener   KSPType        ksp_type;
1152eb910715SAlp Dener 
1153eb910715SAlp Dener   PetscFunctionBegin;
11544f4fdda4SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
11558d5ead36SAlp 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);
11568d5ead36SAlp 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);
11572f75a4aaSAlp 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);
1158748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1159748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1160748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1161748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1162748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1163748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1164748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1165748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1166748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1167748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1168748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1169748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1170748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1171748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1172748696b1SAlp 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);
1173748696b1SAlp 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);
1174748696b1SAlp 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);
1175748696b1SAlp 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);
1176748696b1SAlp 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);
1177748696b1SAlp 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);
1178748696b1SAlp 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);
1179748696b1SAlp 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);
1180748696b1SAlp 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);
1181748696b1SAlp 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);
1182748696b1SAlp 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);
1183748696b1SAlp 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);
1184748696b1SAlp 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);
1185748696b1SAlp 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);
1186748696b1SAlp 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);
1187748696b1SAlp 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);
1188748696b1SAlp 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);
1189748696b1SAlp 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);
1190748696b1SAlp 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);
1191748696b1SAlp 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);
1192748696b1SAlp 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);
1193748696b1SAlp 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);
1194748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1195748696b1SAlp 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);
1196748696b1SAlp 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);
1197748696b1SAlp 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);
1198748696b1SAlp 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);
1199748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1200748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1201748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1202748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1203748696b1SAlp 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);
1204748696b1SAlp 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);
1205c0f10754SAlp 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);
1206eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
120733c78596SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1208eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1209eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1210e0ed867bSAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
121105de396fSBarry Smith   ierr = PetscStrcmp(ksp_type,KSPNASH,&bnk->is_nash);CHKERRQ(ierr);
121205de396fSBarry Smith   ierr = PetscStrcmp(ksp_type,KSPSTCG,&bnk->is_stcg);CHKERRQ(ierr);
121305de396fSBarry Smith   ierr = PetscStrcmp(ksp_type,KSPGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1214eb910715SAlp Dener   PetscFunctionReturn(0);
1215eb910715SAlp Dener }
1216eb910715SAlp Dener 
1217eb910715SAlp Dener /*------------------------------------------------------------*/
1218df278d8fSAlp Dener 
1219e0ed867bSAlp Dener PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1220eb910715SAlp Dener {
1221eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1222eb910715SAlp Dener   PetscInt       nrejects;
1223eb910715SAlp Dener   PetscBool      isascii;
1224eb910715SAlp Dener   PetscErrorCode ierr;
1225eb910715SAlp Dener 
1226eb910715SAlp Dener   PetscFunctionBegin;
1227eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1228eb910715SAlp Dener   if (isascii) {
1229eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1230b9ac7092SAlp Dener     if (bnk->M) {
1231cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(bnk->M,&nrejects);CHKERRQ(ierr);
1232b9ac7092SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %D\n",nrejects);CHKERRQ(ierr);
1233eb910715SAlp Dener     }
1234e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1235eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1236b9ac7092SAlp Dener     if (bnk->M) {
1237eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1238b9ac7092SAlp Dener     }
1239eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1240eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1241eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1242eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1243eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1244eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1245eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1246eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1247eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1248eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1249eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1250eb910715SAlp Dener   }
1251eb910715SAlp Dener   PetscFunctionReturn(0);
1252eb910715SAlp Dener }
1253eb910715SAlp Dener 
1254eb910715SAlp Dener /* ---------------------------------------------------------- */
1255df278d8fSAlp Dener 
1256eb910715SAlp Dener /*MC
1257eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
125866ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1259eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1260eb910715SAlp Dener               Hk dk = -gk
12612b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12622b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1263eb910715SAlp Dener 
1264eb910715SAlp Dener     Options Database Keys:
1265e0ed867bSAlp Dener + -max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1266e0ed867bSAlp Dener . -init_type - trust radius initialization method ("constant", "direction", "interpolation")
1267e0ed867bSAlp Dener . -update_type - trust radius update method ("step", "direction", "interpolation")
1268e0ed867bSAlp Dener . -as_type - active-set estimation method ("none", "bertsekas")
1269e0ed867bSAlp Dener . -as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
1270e0ed867bSAlp Dener . -as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
1271e0ed867bSAlp Dener . -sval - (developer) Hessian perturbation starting value
1272e0ed867bSAlp Dener . -imin - (developer) minimum initial Hessian perturbation
1273e0ed867bSAlp Dener . -imax - (developer) maximum initial Hessian perturbation
1274e0ed867bSAlp Dener . -pmin - (developer) minimum Hessian perturbation
1275e0ed867bSAlp Dener . -pmax - (developer) aximum Hessian perturbation
1276e0ed867bSAlp Dener . -pgfac - (developer) Hessian perturbation growth factor
1277e0ed867bSAlp Dener . -psfac - (developer) Hessian perturbation shrink factor
1278e0ed867bSAlp Dener . -imfac - (developer) initial merit factor for Hessian perturbation
1279e0ed867bSAlp Dener . -pmgfac - (developer) merit growth factor for Hessian perturbation
1280e0ed867bSAlp Dener . -pmsfac - (developer) merit shrink factor for Hessian perturbation
1281e0ed867bSAlp Dener . -eta1 - (developer) threshold for rejecting step (-update_type reduction)
1282e0ed867bSAlp Dener . -eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
1283e0ed867bSAlp Dener . -eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
1284e0ed867bSAlp Dener . -eta4 - (developer) threshold for accepting good step (-update_type reduction)
1285e0ed867bSAlp Dener . -alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
1286e0ed867bSAlp Dener . -alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
1287e0ed867bSAlp Dener . -alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
1288e0ed867bSAlp Dener . -alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
1289e0ed867bSAlp Dener . -alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
1290e0ed867bSAlp Dener . -epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
1291e0ed867bSAlp Dener . -mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
1292e0ed867bSAlp Dener . -mu2 - (developer) threshold for accepting good step (-update_type interpolation)
1293e0ed867bSAlp Dener . -gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
1294e0ed867bSAlp Dener . -gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
1295e0ed867bSAlp Dener . -gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
1296e0ed867bSAlp Dener . -gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
1297e0ed867bSAlp Dener . -theta - (developer) trust region interpolation factor (-update_type interpolation)
1298e0ed867bSAlp Dener . -nu1 - (developer) threshold for small line-search step length (-update_type step)
1299e0ed867bSAlp Dener . -nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
1300e0ed867bSAlp Dener . -nu3 - (developer) threshold for large line-search step length (-update_type step)
1301e0ed867bSAlp Dener . -nu4 - (developer) threshold for very large line-search step length (-update_type step)
1302e0ed867bSAlp Dener . -omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
1303e0ed867bSAlp Dener . -omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
1304e0ed867bSAlp Dener . -omega3 - (developer) radius factor for decent line-search step length (-update_type step)
1305e0ed867bSAlp Dener . -omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
1306e0ed867bSAlp Dener . -omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
1307e0ed867bSAlp Dener . -mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
1308e0ed867bSAlp Dener . -mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
1309e0ed867bSAlp Dener . -gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
1310e0ed867bSAlp Dener . -gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
1311e0ed867bSAlp Dener . -gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
1312e0ed867bSAlp Dener . -gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
1313e0ed867bSAlp Dener - -theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1314eb910715SAlp Dener 
1315eb910715SAlp Dener   Level: beginner
1316eb910715SAlp Dener M*/
1317eb910715SAlp Dener 
1318eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1319eb910715SAlp Dener {
1320eb910715SAlp Dener   TAO_BNK        *bnk;
1321eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1322eb910715SAlp Dener   PetscErrorCode ierr;
1323b9ac7092SAlp Dener   PC             pc;
1324eb910715SAlp Dener 
1325eb910715SAlp Dener   PetscFunctionBegin;
1326eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1327eb910715SAlp Dener 
1328eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1329eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1330eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1331eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1332eb910715SAlp Dener 
1333eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1334eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1335eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1336eb910715SAlp Dener 
1337eb910715SAlp Dener   tao->data = (void*)bnk;
1338eb910715SAlp Dener 
133966ed3702SAlp Dener   /*  Hessian shifting parameters */
1340e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1341e0ed867bSAlp Dener   bnk->computestep = TaoBNKComputeStep;
1342e0ed867bSAlp Dener 
1343eb910715SAlp Dener   bnk->sval   = 0.0;
1344eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1345eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1346eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1347eb910715SAlp Dener 
1348eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1349eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1350eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1351eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1352eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1353eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1354eb910715SAlp Dener 
1355eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1356eb910715SAlp Dener   bnk->nu1 = 0.25;
1357eb910715SAlp Dener   bnk->nu2 = 0.50;
1358eb910715SAlp Dener   bnk->nu3 = 1.00;
1359eb910715SAlp Dener   bnk->nu4 = 1.25;
1360eb910715SAlp Dener 
1361eb910715SAlp Dener   bnk->omega1 = 0.25;
1362eb910715SAlp Dener   bnk->omega2 = 0.50;
1363eb910715SAlp Dener   bnk->omega3 = 1.00;
1364eb910715SAlp Dener   bnk->omega4 = 2.00;
1365eb910715SAlp Dener   bnk->omega5 = 4.00;
1366eb910715SAlp Dener 
1367eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1368eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1369eb910715SAlp Dener   bnk->eta2 = 0.25;
1370eb910715SAlp Dener   bnk->eta3 = 0.50;
1371eb910715SAlp Dener   bnk->eta4 = 0.90;
1372eb910715SAlp Dener 
1373eb910715SAlp Dener   bnk->alpha1 = 0.25;
1374eb910715SAlp Dener   bnk->alpha2 = 0.50;
1375eb910715SAlp Dener   bnk->alpha3 = 1.00;
1376eb910715SAlp Dener   bnk->alpha4 = 2.00;
1377eb910715SAlp Dener   bnk->alpha5 = 4.00;
1378eb910715SAlp Dener 
1379eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1380eb910715SAlp Dener   bnk->mu1 = 0.10;
1381eb910715SAlp Dener   bnk->mu2 = 0.50;
1382eb910715SAlp Dener 
1383eb910715SAlp Dener   bnk->gamma1 = 0.25;
1384eb910715SAlp Dener   bnk->gamma2 = 0.50;
1385eb910715SAlp Dener   bnk->gamma3 = 2.00;
1386eb910715SAlp Dener   bnk->gamma4 = 4.00;
1387eb910715SAlp Dener 
1388eb910715SAlp Dener   bnk->theta = 0.05;
1389eb910715SAlp Dener 
1390eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1391eb910715SAlp Dener   bnk->mu1_i = 0.35;
1392eb910715SAlp Dener   bnk->mu2_i = 0.50;
1393eb910715SAlp Dener 
1394eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1395eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1396eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1397eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1398eb910715SAlp Dener 
1399eb910715SAlp Dener   bnk->theta_i = 0.25;
1400eb910715SAlp Dener 
1401eb910715SAlp Dener   /*  Remaining parameters */
1402c0f10754SAlp Dener   bnk->max_cg_its = 0;
1403eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1404eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1405770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14060a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14070a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
140862675beeSAlp Dener   bnk->dmin = 1.0e-6;
140962675beeSAlp Dener   bnk->dmax = 1.0e6;
1410eb910715SAlp Dener 
1411b9ac7092SAlp Dener   bnk->M               = 0;
1412b9ac7092SAlp Dener   bnk->bfgs_pre        = 0;
1413eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
14147b1c7716SAlp Dener   bnk->update_type     = BNK_UPDATE_REDUCTION;
14152f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1416eb910715SAlp Dener 
1417e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1418e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1419e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1420e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1421e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1422e031d6f5SAlp Dener 
1423c0f10754SAlp Dener   /* Create the line search */
1424eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1425eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1426e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1427eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1428eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1429eb910715SAlp Dener 
1430eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1431eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1432eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1433e0ed867bSAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,"tao_bnk_");CHKERRQ(ierr);
143405de396fSBarry Smith   ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr);
1435f5a7d1c1SBarry Smith   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
1436b9ac7092SAlp Dener   ierr = PCSetType(pc, PCLMVM);CHKERRQ(ierr);
1437eb910715SAlp Dener   PetscFunctionReturn(0);
1438eb910715SAlp Dener }
1439