xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision c5e9d94cca63fa501f774db38e8fb142ef7f5835)
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);
48691b26d3SBarry 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 {
123*c5e9d94cSAlp Dener             ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
124*c5e9d94cSAlp Dener             bnk->H_inactive = tao->hessian;
12528017e9fSAlp Dener           }
126c0f10754SAlp Dener           *needH = PETSC_FALSE;
127eb910715SAlp Dener         }
128eb910715SAlp Dener 
129eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
13062602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
13162602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
13262602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
1333b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
13489da521bSAlp Dener           /* Compute the step we actually accepted */
135eb910715SAlp Dener           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
13662602cfbSAlp Dener           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
13762602cfbSAlp Dener           /* Compute the objective at the trial */
13862602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
139691b26d3SBarry Smith           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
14062602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
141eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
142eb910715SAlp Dener             tau = bnk->gamma1_i;
143eb910715SAlp Dener           } else {
1440a4511e9SAlp Dener             if (ftrial < f_min) {
1450a4511e9SAlp Dener               f_min = ftrial;
146eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
147eb910715SAlp Dener             }
14808752603SAlp Dener 
149770b7498SAlp Dener             /* Compute the predicted and actual reduction */
15089da521bSAlp Dener             if (bnk->active_idx) {
15108752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
15208752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1532ab2a32cSAlp Dener             } else {
15408752603SAlp Dener               bnk->X_inactive = bnk->W;
15508752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1562ab2a32cSAlp Dener             }
15708752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
15808752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
15989da521bSAlp Dener             if (bnk->active_idx) {
16008752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
16108752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1622ab2a32cSAlp Dener             }
163eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
164eb910715SAlp Dener             actred = bnk->f - ftrial;
1653105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
166eb910715SAlp Dener               kappa = 1.0;
1673105154fSTodd Munson             } else {
168eb910715SAlp Dener               kappa = actred / prered;
169eb910715SAlp Dener             }
170eb910715SAlp Dener 
171eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
172eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
173eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
174eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
175eb910715SAlp Dener 
17618cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
177eb910715SAlp Dener               /*  Great agreement */
178eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
179eb910715SAlp Dener 
180eb910715SAlp Dener               if (tau_max < 1.0) {
181eb910715SAlp Dener                 tau = bnk->gamma3_i;
1823105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
183eb910715SAlp Dener                 tau = bnk->gamma4_i;
1843105154fSTodd Munson               } else {
185eb910715SAlp Dener                 tau = tau_max;
186eb910715SAlp Dener               }
18718cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
188eb910715SAlp Dener               /*  Good agreement */
189eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
190eb910715SAlp Dener 
191eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
192eb910715SAlp Dener                 tau = bnk->gamma2_i;
193eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
194eb910715SAlp Dener                 tau = bnk->gamma3_i;
195eb910715SAlp Dener               } else {
196eb910715SAlp Dener                 tau = tau_max;
197eb910715SAlp Dener               }
1988f8a4e06SAlp Dener             } else {
199eb910715SAlp Dener               /*  Not good agreement */
200eb910715SAlp Dener               if (tau_min > 1.0) {
201eb910715SAlp Dener                 tau = bnk->gamma2_i;
202eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
203eb910715SAlp Dener                 tau = bnk->gamma1_i;
204eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
205eb910715SAlp Dener                 tau = bnk->gamma1_i;
2063105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
207eb910715SAlp Dener                 tau = tau_1;
2083105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
209eb910715SAlp Dener                 tau = tau_2;
210eb910715SAlp Dener               } else {
211eb910715SAlp Dener                 tau = tau_max;
212eb910715SAlp Dener               }
213eb910715SAlp Dener             }
214eb910715SAlp Dener           }
215eb910715SAlp Dener           tao->trust = tau * tao->trust;
216eb910715SAlp Dener         }
217eb910715SAlp Dener 
2180a4511e9SAlp Dener         if (f_min < bnk->f) {
219937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2200a4511e9SAlp Dener           bnk->f = f_min;
221937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
222eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
2233b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
224937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
225937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
22609164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
22761be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
22861be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
22961be54a6SAlp Dener           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
230937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
231f5766c09SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
232c0f10754SAlp Dener           *needH = PETSC_TRUE;
233937a31a1SAlp Dener           /* Test the new step for convergence */
23489da521bSAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
23589da521bSAlp Dener           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
236691b26d3SBarry Smith           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
237e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
238e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
239eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
240eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
241937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
242414d97d3SAlp Dener           ierr = TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
243eb910715SAlp Dener         }
244eb910715SAlp Dener       }
245eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
246e031d6f5SAlp Dener 
247e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
248e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
249e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
250eb910715SAlp Dener       break;
251eb910715SAlp Dener 
252eb910715SAlp Dener     default:
253eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
254eb910715SAlp Dener       tao->trust = 0.0;
255eb910715SAlp Dener       break;
256eb910715SAlp Dener     }
257eb910715SAlp Dener   }
258eb910715SAlp Dener   PetscFunctionReturn(0);
259eb910715SAlp Dener }
260eb910715SAlp Dener 
261df278d8fSAlp Dener /*------------------------------------------------------------*/
262df278d8fSAlp Dener 
263e0ed867bSAlp Dener /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */
26462675beeSAlp Dener 
26562675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
26662675beeSAlp Dener {
26762675beeSAlp Dener   PetscErrorCode               ierr;
26862675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
26962675beeSAlp Dener 
27062675beeSAlp Dener   PetscFunctionBegin;
27162675beeSAlp Dener   /* Compute the Hessian */
27262675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
27362675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
274b9ac7092SAlp Dener   if (bnk->M) {
27562675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
27662675beeSAlp Dener   }
277e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
278f5766c09SAlp Dener   if (bnk->Hpre_inactive) {
279f5766c09SAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
280f5766c09SAlp Dener   }
281f5766c09SAlp Dener   if (bnk->H_inactive) {
282e0ed867bSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
283f5766c09SAlp Dener   }
284f5766c09SAlp Dener   if (bnk->active_idx) {
285e0ed867bSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
286e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
287f5766c09SAlp Dener       ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
288e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
289e0ed867bSAlp Dener     } else {
290e0ed867bSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
291e0ed867bSAlp Dener     }
292e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
293e0ed867bSAlp Dener       ierr = PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx);CHKERRQ(ierr);
294e0ed867bSAlp Dener     }
295e0ed867bSAlp Dener   } else {
296*c5e9d94cSAlp Dener     ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
297*c5e9d94cSAlp Dener     bnk->H_inactive = tao->hessian;
298e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
299f5766c09SAlp Dener       ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
300e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
301e0ed867bSAlp Dener     } else {
302*c5e9d94cSAlp Dener       ierr = PetscObjectReference((PetscObject)tao->hessian_pre);
303*c5e9d94cSAlp Dener       bnk->Hpre_inactive = tao->hessian_pre;
304e0ed867bSAlp Dener     }
305e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
306e0ed867bSAlp Dener       ierr = PCLMVMClearIS(bnk->bfgs_pre);CHKERRQ(ierr);
307e0ed867bSAlp Dener     }
308e0ed867bSAlp Dener   }
30962675beeSAlp Dener   PetscFunctionReturn(0);
31062675beeSAlp Dener }
31162675beeSAlp Dener 
31262675beeSAlp Dener /*------------------------------------------------------------*/
31362675beeSAlp Dener 
3142f75a4aaSAlp Dener /* Routine for estimating the active set */
3152f75a4aaSAlp Dener 
31608752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3172f75a4aaSAlp Dener {
3182f75a4aaSAlp Dener   PetscErrorCode               ierr;
3192f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
320c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3212f75a4aaSAlp Dener 
3222f75a4aaSAlp Dener   PetscFunctionBegin;
32308752603SAlp Dener   switch (asType) {
3242f75a4aaSAlp Dener   case BNK_AS_NONE:
3252f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3262f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3272f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3282f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3292f75a4aaSAlp Dener     break;
3302f75a4aaSAlp Dener 
3312f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3322f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
333b9ac7092SAlp Dener     if (bnk->M) {
3342f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3359515a401SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3362f75a4aaSAlp Dener     } else {
337fc5ca067SStefano Zampini       hessComputed = diagExists = PETSC_FALSE;
338f5766c09SAlp Dener       if (tao->hessian) {
33961be54a6SAlp Dener         ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
340f5766c09SAlp Dener       }
341fc5ca067SStefano Zampini       if (hessComputed) {
342fc5ca067SStefano Zampini         ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
343fc5ca067SStefano Zampini       }
344fc5ca067SStefano Zampini       if (diagExists) {
3459b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3462f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3479b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3489b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
349c3b366b1Sprj-         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);
3502f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
35161be54a6SAlp Dener       } else {
352c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
35361be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
35461be54a6SAlp Dener       }
3552f75a4aaSAlp Dener     }
3562f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
35789da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
35889da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
359c4b75bccSAlp Dener     break;
3602f75a4aaSAlp Dener 
3612f75a4aaSAlp Dener   default:
3622f75a4aaSAlp Dener     break;
3632f75a4aaSAlp Dener   }
3642f75a4aaSAlp Dener   PetscFunctionReturn(0);
3652f75a4aaSAlp Dener }
3662f75a4aaSAlp Dener 
3672f75a4aaSAlp Dener /*------------------------------------------------------------*/
3682f75a4aaSAlp Dener 
3692f75a4aaSAlp Dener /* Routine for bounding the step direction */
3702f75a4aaSAlp Dener 
371a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
3722f75a4aaSAlp Dener {
3732f75a4aaSAlp Dener   PetscErrorCode               ierr;
3742f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3752f75a4aaSAlp Dener 
3762f75a4aaSAlp Dener   PetscFunctionBegin;
377a1318120SAlp Dener   switch (asType) {
3782f75a4aaSAlp Dener   case BNK_AS_NONE:
379c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
3802f75a4aaSAlp Dener     break;
3812f75a4aaSAlp Dener 
3822f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
383c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
3842f75a4aaSAlp Dener     break;
3852f75a4aaSAlp Dener 
3862f75a4aaSAlp Dener   default:
3872f75a4aaSAlp Dener     break;
3882f75a4aaSAlp Dener   }
3892f75a4aaSAlp Dener   PetscFunctionReturn(0);
3902f75a4aaSAlp Dener }
3912f75a4aaSAlp Dener 
392e031d6f5SAlp Dener /*------------------------------------------------------------*/
393e031d6f5SAlp Dener 
394e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
395e031d6f5SAlp Dener    accelerate Newton convergence.
396e031d6f5SAlp Dener 
397e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
398e031d6f5SAlp Dener    for more gradient evaluations.
399e031d6f5SAlp Dener */
400e031d6f5SAlp Dener 
401c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
402c0f10754SAlp Dener {
403c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
404c0f10754SAlp Dener   PetscErrorCode               ierr;
405c0f10754SAlp Dener 
406c0f10754SAlp Dener   PetscFunctionBegin;
407c0f10754SAlp Dener   *terminate = PETSC_FALSE;
408c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
409c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
410c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
411c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
412c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
413c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
414c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
415c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
416c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
417c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
418e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
419c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
420c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
421c0f10754SAlp 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) {
422c0f10754SAlp Dener       *terminate = PETSC_TRUE;
42361be54a6SAlp Dener     } else {
42433c78596SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
425c0f10754SAlp Dener     }
426c0f10754SAlp Dener   }
427c0f10754SAlp Dener   PetscFunctionReturn(0);
428c0f10754SAlp Dener }
429c0f10754SAlp Dener 
4302f75a4aaSAlp Dener /*------------------------------------------------------------*/
4312f75a4aaSAlp Dener 
432c0f10754SAlp Dener /* Routine for computing the Newton step. */
433df278d8fSAlp Dener 
4346b591159SAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
435eb910715SAlp Dener {
436eb910715SAlp Dener   PetscErrorCode               ierr;
437eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
438eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
439eb910715SAlp Dener   PetscInt                     kspits;
440bddd1ffdSAlp Dener   PetscBool                    is_lmvm;
441eb910715SAlp Dener 
442eb910715SAlp Dener   PetscFunctionBegin;
44389da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
44489da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
44589da521bSAlp Dener   if (!bnk->inactive_idx) {
44689da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
447a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
44889da521bSAlp Dener     PetscFunctionReturn(0);
44989da521bSAlp Dener   }
45089da521bSAlp Dener 
45162675beeSAlp Dener   /* Shift the reduced Hessian matrix */
45262675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
453f7bf01afSAlp Dener     ierr = PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm);CHKERRQ(ierr);
454f7bf01afSAlp Dener     if (is_lmvm) {
455f7bf01afSAlp Dener       ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
456f7bf01afSAlp Dener     } else {
45762675beeSAlp Dener       ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
45862675beeSAlp Dener       if (bnk->H_inactive != bnk->Hpre_inactive) {
45962675beeSAlp Dener         ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
46062675beeSAlp Dener       }
46162675beeSAlp Dener     }
462f7bf01afSAlp Dener   }
46362675beeSAlp Dener 
464eb910715SAlp Dener   /* Solve the Newton system of equations */
465937a31a1SAlp Dener   tao->ksp_its = 0;
4662f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
4675e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
46894d5fdc2STristan Konolige   ierr = KSPResetFromOptions(tao->ksp);CHKERRQ(ierr);
46909164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
4705e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
47189da521bSAlp Dener   if (bnk->active_idx) {
4725e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
4735e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
4745e9b73cbSAlp Dener   } else {
4755e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4765e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
47728017e9fSAlp Dener   }
478eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
479fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4805e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
481eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
482eb910715SAlp Dener     tao->ksp_its+=kspits;
483eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
484080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
485eb910715SAlp Dener 
486eb910715SAlp Dener     if (0.0 == tao->trust) {
487eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
488080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
489080d2917SAlp Dener         tao->trust = bnk->dnorm;
490eb910715SAlp Dener 
491eb910715SAlp Dener         /* Modify the radius if it is too large or small */
492eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
493eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
494eb910715SAlp Dener       } else {
495eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
496eb910715SAlp Dener            the trust-region subproblem to get a direction */
497eb910715SAlp Dener         tao->trust = tao->trust0;
498eb910715SAlp Dener 
499eb910715SAlp Dener         /* Modify the radius if it is too large or small */
500eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
501eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
502eb910715SAlp Dener 
503fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5045e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
505eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
506eb910715SAlp Dener         tao->ksp_its+=kspits;
507eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
508080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
509eb910715SAlp Dener 
510691b26d3SBarry Smith         if (bnk->dnorm == 0.0) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero");
511eb910715SAlp Dener       }
512eb910715SAlp Dener     }
513eb910715SAlp Dener   } else {
5145e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
515eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
516eb910715SAlp Dener     tao->ksp_its += kspits;
517eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
518eb910715SAlp Dener   }
5195e9b73cbSAlp Dener   /* Restore sub vectors back */
52089da521bSAlp Dener   if (bnk->active_idx) {
5215e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5225e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5235e9b73cbSAlp Dener   }
524770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
525fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
526a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
527770b7498SAlp Dener 
528770b7498SAlp Dener   /* Record convergence reasons */
529e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
530e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
531770b7498SAlp Dener     ++bnk->ksp_atol;
532e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
533770b7498SAlp Dener     ++bnk->ksp_rtol;
534e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
535770b7498SAlp Dener     ++bnk->ksp_ctol;
536e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
537770b7498SAlp Dener     ++bnk->ksp_negc;
538e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
539770b7498SAlp Dener     ++bnk->ksp_dtol;
540e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
541770b7498SAlp Dener     ++bnk->ksp_iter;
542770b7498SAlp Dener   } else {
543770b7498SAlp Dener     ++bnk->ksp_othr;
544770b7498SAlp Dener   }
545fed79b8eSAlp Dener 
546fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
547b9ac7092SAlp Dener   if (bnk->M) {
548cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
549b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
550fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
551cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
55209164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
553eb910715SAlp Dener     }
554fed79b8eSAlp Dener   }
5556b591159SAlp Dener   *step_type = BNK_NEWTON;
556e465cd6fSAlp Dener   PetscFunctionReturn(0);
557e465cd6fSAlp Dener }
558eb910715SAlp Dener 
55962675beeSAlp Dener /*------------------------------------------------------------*/
56062675beeSAlp Dener 
5615e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5625e9b73cbSAlp Dener 
5635e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
5645e9b73cbSAlp Dener {
5655e9b73cbSAlp Dener   PetscErrorCode               ierr;
5665e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
5675e9b73cbSAlp Dener 
5685e9b73cbSAlp Dener   PetscFunctionBegin;
5695e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
57089da521bSAlp Dener   if (bnk->active_idx) {
5715e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5725e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5735e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5745e9b73cbSAlp Dener   } else {
5755e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
5765e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5775e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
5785e9b73cbSAlp Dener   }
5795e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5805e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
5815e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
58233c78596SAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);CHKERRQ(ierr);
5835e9b73cbSAlp Dener   /* Restore the sub vectors */
58489da521bSAlp Dener   if (bnk->active_idx) {
5855e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5865e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5875e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5885e9b73cbSAlp Dener   }
5895e9b73cbSAlp Dener   PetscFunctionReturn(0);
5905e9b73cbSAlp Dener }
5915e9b73cbSAlp Dener 
5925e9b73cbSAlp Dener /*------------------------------------------------------------*/
5935e9b73cbSAlp Dener 
59462675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
59562675beeSAlp Dener 
59662675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
59762675beeSAlp Dener    in the event that the Newton step fails the test.
59862675beeSAlp Dener */
59962675beeSAlp Dener 
600e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
601e465cd6fSAlp Dener {
602e465cd6fSAlp Dener   PetscErrorCode               ierr;
603e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
604e465cd6fSAlp Dener 
605b2d8c577SAlp Dener   PetscReal                    gdx, e_min;
606e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
607e465cd6fSAlp Dener 
608e465cd6fSAlp Dener   PetscFunctionBegin;
6096b591159SAlp Dener   switch (*stepType) {
6106b591159SAlp Dener   case BNK_NEWTON:
611080d2917SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
612eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
613eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
614eb910715SAlp Dener         Update the perturbation for next time */
615eb910715SAlp Dener       if (bnk->pert <= 0.0) {
616eb910715SAlp Dener         /* Initialize the perturbation */
617eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
618eb910715SAlp Dener         if (bnk->is_gltr) {
61905de396fSBarry Smith           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
620eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
621eb910715SAlp Dener         }
622eb910715SAlp Dener       } else {
623eb910715SAlp Dener         /* Increase the perturbation */
624eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
625eb910715SAlp Dener       }
626eb910715SAlp Dener 
6270ad3a497SAlp Dener       if (!bnk->M) {
628eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
629eb910715SAlp Dener           Must use gradient direction in this case */
630080d2917SAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
631eb910715SAlp Dener         *stepType = BNK_GRADIENT;
632eb910715SAlp Dener       } else {
633eb910715SAlp Dener         /* Attempt to use the BFGS direction */
6349515a401SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
635eb910715SAlp Dener 
6368d5ead36SAlp Dener         /* Check for success (descent direction)
6378d5ead36SAlp Dener           NOTE: Negative gdx here means not a descent direction because
6388d5ead36SAlp Dener           the fall-back step is missing a negative sign. */
639080d2917SAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
6403105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
641eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
642eb910715SAlp Dener             We can assert bfgsUpdates > 1 in this case because
643eb910715SAlp Dener             the first solve produces the scaled gradient direction,
644eb910715SAlp Dener             which is guaranteed to be descent */
645eb910715SAlp Dener 
646eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
647cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
64809164190SAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
6499515a401SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
650eb910715SAlp Dener 
651eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
652eb910715SAlp Dener         } else {
653cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
654eb910715SAlp Dener           if (1 == bfgsUpdates) {
655eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
656eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
657eb910715SAlp Dener           } else {
658eb910715SAlp Dener             *stepType = BNK_BFGS;
659eb910715SAlp Dener           }
660eb910715SAlp Dener         }
661eb910715SAlp Dener       }
6628d5ead36SAlp Dener       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6638d5ead36SAlp Dener       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
664a1318120SAlp Dener       ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
665eb910715SAlp Dener     } else {
666eb910715SAlp Dener       /* Computed Newton step is descent */
667eb910715SAlp Dener       switch (ksp_reason) {
668eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
669eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
670eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
671eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
672eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
673eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
674eb910715SAlp Dener         if (bnk->pert <= 0.0) {
675eb910715SAlp Dener           /* Initialize the perturbation */
676eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
677eb910715SAlp Dener           if (bnk->is_gltr) {
67805de396fSBarry Smith             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
679eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
680eb910715SAlp Dener           }
681eb910715SAlp Dener         } else {
682eb910715SAlp Dener           /* Increase the perturbation */
683eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
684eb910715SAlp Dener         }
685eb910715SAlp Dener         break;
686eb910715SAlp Dener 
687eb910715SAlp Dener       default:
688eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
689eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
690eb910715SAlp Dener         if (bnk->pert < bnk->pmin) {
691eb910715SAlp Dener           bnk->pert = 0.0;
692eb910715SAlp Dener         }
693eb910715SAlp Dener         break;
694eb910715SAlp Dener       }
695fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
696eb910715SAlp Dener     }
6976b591159SAlp Dener     break;
6986b591159SAlp Dener 
6996b591159SAlp Dener   case BNK_BFGS:
7006b591159SAlp Dener     /* Check for success (descent direction) */
7016b591159SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
7026b591159SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
7036b591159SAlp Dener       /* Step is not descent or solve was not successful
7046b591159SAlp Dener          Use steepest descent direction (scaled) */
7056b591159SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
7066b591159SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
7079515a401SAlp Dener       ierr = MatSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
7086b591159SAlp Dener       ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
7096b591159SAlp Dener       ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
7106b591159SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
7116b591159SAlp Dener     } else {
7126b591159SAlp Dener       *stepType = BNK_BFGS;
7136b591159SAlp Dener     }
7146b591159SAlp Dener     break;
7156b591159SAlp Dener 
7166b591159SAlp Dener   case BNK_SCALED_GRADIENT:
7176b591159SAlp Dener     break;
7186b591159SAlp Dener 
7196b591159SAlp Dener   default:
7206b591159SAlp Dener     break;
7216b591159SAlp Dener   }
7226b591159SAlp Dener 
723eb910715SAlp Dener   PetscFunctionReturn(0);
724eb910715SAlp Dener }
725eb910715SAlp Dener 
726df278d8fSAlp Dener /*------------------------------------------------------------*/
727df278d8fSAlp Dener 
728df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
729df278d8fSAlp Dener 
730df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
731df278d8fSAlp Dener   Newton step does not produce a valid step length.
732df278d8fSAlp Dener */
733df278d8fSAlp Dener 
734937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
735c14b763aSAlp Dener {
736c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
737c14b763aSAlp Dener   PetscErrorCode ierr;
738c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
739c14b763aSAlp Dener 
740b2d8c577SAlp Dener   PetscReal      e_min, gdx;
741c14b763aSAlp Dener   PetscInt       bfgsUpdates;
742c14b763aSAlp Dener 
743c14b763aSAlp Dener   PetscFunctionBegin;
744c14b763aSAlp Dener   /* Perform the linesearch */
745c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
746c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
747c14b763aSAlp Dener 
748b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
749c14b763aSAlp Dener     /* Linesearch failed, revert solution */
750c14b763aSAlp Dener     bnk->f = bnk->fold;
751c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
752c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
753c14b763aSAlp Dener 
754937a31a1SAlp Dener     switch(*stepType) {
755c14b763aSAlp Dener     case BNK_NEWTON:
7568d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
757c14b763aSAlp Dener          Update the perturbation for next time */
758c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
759c14b763aSAlp Dener         /* Initialize the perturbation */
760c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
761c14b763aSAlp Dener         if (bnk->is_gltr) {
76205de396fSBarry Smith           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
763c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
764c14b763aSAlp Dener         }
765c14b763aSAlp Dener       } else {
766c14b763aSAlp Dener         /* Increase the perturbation */
767c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
768c14b763aSAlp Dener       }
769c14b763aSAlp Dener 
7700ad3a497SAlp Dener       if (!bnk->M) {
771c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
772c14b763aSAlp Dener            Must use gradient direction in this case */
773937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
774937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
775c14b763aSAlp Dener       } else {
776c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
7779515a401SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
7788d5ead36SAlp Dener         /* Check for success (descent direction)
7798d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
780c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7813105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
782c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
783c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
784c14b763aSAlp Dener              Use steepest descent direction (scaled) */
785cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
786c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
7879515a401SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
788c14b763aSAlp Dener 
789c14b763aSAlp Dener           bfgsUpdates = 1;
790937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
791c14b763aSAlp Dener         } else {
792cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
793c14b763aSAlp Dener           if (1 == bfgsUpdates) {
794c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
795937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
796c14b763aSAlp Dener           } else {
797937a31a1SAlp Dener             *stepType = BNK_BFGS;
798c14b763aSAlp Dener           }
799c14b763aSAlp Dener         }
800c14b763aSAlp Dener       }
801c14b763aSAlp Dener       break;
802c14b763aSAlp Dener 
803c14b763aSAlp Dener     case BNK_BFGS:
804c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
805c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
806c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
807cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
808c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
8099515a401SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
810c14b763aSAlp Dener 
811c14b763aSAlp Dener       bfgsUpdates = 1;
812937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
813c14b763aSAlp Dener       break;
814c14b763aSAlp Dener     }
8158d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8168d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
817a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
818c14b763aSAlp Dener 
8198d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
820c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
821c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
822c14b763aSAlp Dener   }
823c14b763aSAlp Dener   *reason = ls_reason;
824c14b763aSAlp Dener   PetscFunctionReturn(0);
825c14b763aSAlp Dener }
826c14b763aSAlp Dener 
827df278d8fSAlp Dener /*------------------------------------------------------------*/
828df278d8fSAlp Dener 
829df278d8fSAlp Dener /* Routine for updating the trust radius.
830df278d8fSAlp Dener 
831df278d8fSAlp Dener   Function features three different update methods:
832df278d8fSAlp Dener   1) Line-search step length based
833df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
834df278d8fSAlp Dener   3) Interpolation
835df278d8fSAlp Dener */
836df278d8fSAlp Dener 
83728017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
838080d2917SAlp Dener {
839080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
840080d2917SAlp Dener   PetscErrorCode ierr;
841080d2917SAlp Dener 
842b1c2d0e3SAlp Dener   PetscReal      step, kappa;
843080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
844080d2917SAlp Dener 
845080d2917SAlp Dener   PetscFunctionBegin;
846080d2917SAlp Dener   /* Update trust region radius */
847080d2917SAlp Dener   *accept = PETSC_FALSE;
84828017e9fSAlp Dener   switch(updateType) {
849080d2917SAlp Dener   case BNK_UPDATE_STEP:
850c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
851080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
852080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
853080d2917SAlp Dener       if (step < bnk->nu1) {
854080d2917SAlp Dener         /* Very bad step taken; reduce radius */
855080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
856080d2917SAlp Dener       } else if (step < bnk->nu2) {
857080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
858080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
859080d2917SAlp Dener       } else if (step < bnk->nu3) {
860080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
861080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
862080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
863080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
864080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
865080d2917SAlp Dener         }
866080d2917SAlp Dener       } else if (step < bnk->nu4) {
867080d2917SAlp Dener         /*  Full step taken; increase the radius */
868080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
869080d2917SAlp Dener       } else {
870080d2917SAlp Dener         /*  More than full step taken; increase the radius */
871080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
872080d2917SAlp Dener       }
873080d2917SAlp Dener     } else {
874080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
875080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
876080d2917SAlp Dener     }
877080d2917SAlp Dener     break;
878080d2917SAlp Dener 
879080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
880080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
881e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
882fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
883fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
884fed79b8eSAlp Dener            be rejected! */
885080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8863105154fSTodd Munson       } else {
887b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
888080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
889080d2917SAlp Dener         } else {
8903105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
891080d2917SAlp Dener             kappa = 1.0;
8923105154fSTodd Munson           } else {
893080d2917SAlp Dener             kappa = actred / prered;
894080d2917SAlp Dener           }
895fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
896080d2917SAlp Dener           if (kappa < bnk->eta1) {
897fed79b8eSAlp Dener             /* Reject the step */
898080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8993105154fSTodd Munson           } else {
900fed79b8eSAlp Dener             /* Accept the step */
901c133c014SAlp Dener             *accept = PETSC_TRUE;
902c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
9038d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
904080d2917SAlp Dener               if (kappa < bnk->eta2) {
905080d2917SAlp Dener                 /* Marginal bad step */
906c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
9073105154fSTodd Munson               } else if (kappa < bnk->eta3) {
908fed79b8eSAlp Dener                 /* Reasonable step */
909fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
9103105154fSTodd Munson               } else if (kappa < bnk->eta4) {
911080d2917SAlp Dener                 /* Good step */
912c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
9133105154fSTodd Munson               } else {
914080d2917SAlp Dener                 /* Very good step */
915c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
916080d2917SAlp Dener               }
917c133c014SAlp Dener             }
918080d2917SAlp Dener           }
919080d2917SAlp Dener         }
920080d2917SAlp Dener       }
921080d2917SAlp Dener     } else {
922080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
923080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
924080d2917SAlp Dener     }
925080d2917SAlp Dener     break;
926080d2917SAlp Dener 
927080d2917SAlp Dener   default:
928080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
929b1c2d0e3SAlp Dener       if (prered < 0.0) {
930080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
931080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
932080d2917SAlp Dener         /*  be rejected! */
933080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
934080d2917SAlp Dener       } else {
935b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
936080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
937080d2917SAlp Dener         } else {
938080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
939080d2917SAlp Dener             kappa = 1.0;
940080d2917SAlp Dener           } else {
941080d2917SAlp Dener             kappa = actred / prered;
942080d2917SAlp Dener           }
943080d2917SAlp Dener 
944080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
945080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
946080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
947080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
948080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
949080d2917SAlp Dener 
950080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
951080d2917SAlp Dener             /*  Great agreement */
952080d2917SAlp Dener             *accept = PETSC_TRUE;
953080d2917SAlp Dener             if (tau_max < 1.0) {
954080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
955080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
956080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
957080d2917SAlp Dener             } else {
958080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
959080d2917SAlp Dener             }
960080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
961080d2917SAlp Dener             /*  Good agreement */
962080d2917SAlp Dener             *accept = PETSC_TRUE;
963080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
964080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
965080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
966080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
967080d2917SAlp Dener             } else if (tau_max < 1.0) {
968080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
969080d2917SAlp Dener             } else {
970080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
971080d2917SAlp Dener             }
972080d2917SAlp Dener           } else {
973080d2917SAlp Dener             /*  Not good agreement */
974080d2917SAlp Dener             if (tau_min > 1.0) {
975080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
976080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
977080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
978080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
979080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
980080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
981080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
982080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
983080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
984080d2917SAlp Dener             } else {
985080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
986080d2917SAlp Dener             }
987080d2917SAlp Dener           }
988080d2917SAlp Dener         }
989080d2917SAlp Dener       }
990080d2917SAlp Dener     } else {
991080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
992080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
993080d2917SAlp Dener     }
99428017e9fSAlp Dener     break;
995080d2917SAlp Dener   }
996c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
997c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
998fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
999080d2917SAlp Dener   PetscFunctionReturn(0);
1000080d2917SAlp Dener }
1001080d2917SAlp Dener 
1002eb910715SAlp Dener /* ---------------------------------------------------------- */
1003df278d8fSAlp Dener 
100462675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
100562675beeSAlp Dener {
100662675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
100762675beeSAlp Dener 
100862675beeSAlp Dener   PetscFunctionBegin;
100962675beeSAlp Dener   switch (stepType) {
101062675beeSAlp Dener   case BNK_NEWTON:
101162675beeSAlp Dener     ++bnk->newt;
101262675beeSAlp Dener     break;
101362675beeSAlp Dener   case BNK_BFGS:
101462675beeSAlp Dener     ++bnk->bfgs;
101562675beeSAlp Dener     break;
101662675beeSAlp Dener   case BNK_SCALED_GRADIENT:
101762675beeSAlp Dener     ++bnk->sgrad;
101862675beeSAlp Dener     break;
101962675beeSAlp Dener   case BNK_GRADIENT:
102062675beeSAlp Dener     ++bnk->grad;
102162675beeSAlp Dener     break;
102262675beeSAlp Dener   default:
102362675beeSAlp Dener     break;
102462675beeSAlp Dener   }
102562675beeSAlp Dener   PetscFunctionReturn(0);
102662675beeSAlp Dener }
102762675beeSAlp Dener 
102862675beeSAlp Dener /* ---------------------------------------------------------- */
102962675beeSAlp Dener 
10309b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1031eb910715SAlp Dener {
1032eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1033eb910715SAlp Dener   PetscErrorCode ierr;
1034e031d6f5SAlp Dener   PetscInt       i;
10355eb5f4d6SAlp Dener   KSPType        ksp_type;
1036eb910715SAlp Dener 
1037eb910715SAlp Dener   PetscFunctionBegin;
1038c4b75bccSAlp Dener   if (!tao->gradient) {
1039c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1040c4b75bccSAlp Dener   }
1041c4b75bccSAlp Dener   if (!tao->stepdirection) {
1042c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1043c4b75bccSAlp Dener   }
1044c4b75bccSAlp Dener   if (!bnk->W) {
1045c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1046c4b75bccSAlp Dener   }
1047c4b75bccSAlp Dener   if (!bnk->Xold) {
1048c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1049c4b75bccSAlp Dener   }
1050c4b75bccSAlp Dener   if (!bnk->Gold) {
1051c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1052c4b75bccSAlp Dener   }
1053c4b75bccSAlp Dener   if (!bnk->Xwork) {
1054c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1055c4b75bccSAlp Dener   }
1056c4b75bccSAlp Dener   if (!bnk->Gwork) {
1057c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1058c4b75bccSAlp Dener   }
1059c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1060c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1061c4b75bccSAlp Dener   }
1062c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1063c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1064c4b75bccSAlp Dener   }
1065c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1066c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1067c4b75bccSAlp Dener   }
1068c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1069c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1070c4b75bccSAlp Dener   }
1071e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1072c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1073c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
107489da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
107589da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
107689da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1077c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1078c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1079c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1080c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1081c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1082c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1083c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1084c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1085c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1086c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1087c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1088c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1089c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1090c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1091e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1092e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1093e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1094937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1095e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1096e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1097e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1098e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1099c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1100e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1101e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1102e031d6f5SAlp Dener     }
1103e031d6f5SAlp Dener   }
11045eb5f4d6SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
11055eb5f4d6SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPNASH,&bnk->is_nash);CHKERRQ(ierr);
11065eb5f4d6SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPSTCG,&bnk->is_stcg);CHKERRQ(ierr);
11075eb5f4d6SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPGLTR,&bnk->is_gltr);CHKERRQ(ierr);
110883c8fe1dSLisandro Dalcin   bnk->X_inactive = NULL;
110983c8fe1dSLisandro Dalcin   bnk->G_inactive = NULL;
111083c8fe1dSLisandro Dalcin   bnk->inactive_work = NULL;
111183c8fe1dSLisandro Dalcin   bnk->active_work = NULL;
111283c8fe1dSLisandro Dalcin   bnk->inactive_idx = NULL;
111383c8fe1dSLisandro Dalcin   bnk->active_idx = NULL;
111483c8fe1dSLisandro Dalcin   bnk->active_lower = NULL;
111583c8fe1dSLisandro Dalcin   bnk->active_upper = NULL;
111683c8fe1dSLisandro Dalcin   bnk->active_fixed = NULL;
111783c8fe1dSLisandro Dalcin   bnk->M = NULL;
111883c8fe1dSLisandro Dalcin   bnk->H_inactive = NULL;
111983c8fe1dSLisandro Dalcin   bnk->Hpre_inactive = NULL;
1120eb910715SAlp Dener   PetscFunctionReturn(0);
1121eb910715SAlp Dener }
1122eb910715SAlp Dener 
1123eb910715SAlp Dener /*------------------------------------------------------------*/
1124df278d8fSAlp Dener 
1125e0ed867bSAlp Dener PetscErrorCode TaoDestroy_BNK(Tao tao)
1126eb910715SAlp Dener {
1127eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1128eb910715SAlp Dener   PetscErrorCode ierr;
1129eb910715SAlp Dener 
1130eb910715SAlp Dener   PetscFunctionBegin;
1131eb910715SAlp Dener   if (tao->setupcalled) {
1132eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1133eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1134eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
113509164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
113609164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
113709164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
113809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
113962675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
114062675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1141c4b75bccSAlp Dener   }
1142ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1143ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1144ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1145ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1146ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
1147c4b75bccSAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1148c4b75bccSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1149ca964c49SAlp Dener   ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1150eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1151eb910715SAlp Dener   PetscFunctionReturn(0);
1152eb910715SAlp Dener }
1153eb910715SAlp Dener 
1154eb910715SAlp Dener /*------------------------------------------------------------*/
1155df278d8fSAlp Dener 
1156e0ed867bSAlp Dener PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1157eb910715SAlp Dener {
1158eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1159eb910715SAlp Dener   PetscErrorCode ierr;
1160eb910715SAlp Dener 
1161eb910715SAlp Dener   PetscFunctionBegin;
11624f4fdda4SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
116383c8fe1dSLisandro Dalcin   ierr = PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL);CHKERRQ(ierr);
116483c8fe1dSLisandro Dalcin   ierr = PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL);CHKERRQ(ierr);
116583c8fe1dSLisandro Dalcin   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL);CHKERRQ(ierr);
1166748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1167748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1168748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1169748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1170748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1171748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1172748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1173748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1174748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1175748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1176748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1177748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1178748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1179748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1180748696b1SAlp 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);
1181748696b1SAlp 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);
1182748696b1SAlp 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);
1183748696b1SAlp 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);
1184748696b1SAlp 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);
1185748696b1SAlp 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);
1186748696b1SAlp 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);
1187748696b1SAlp 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);
1188748696b1SAlp 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);
1189748696b1SAlp 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);
1190748696b1SAlp 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);
1191748696b1SAlp 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);
1192748696b1SAlp 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);
1193748696b1SAlp 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);
1194748696b1SAlp 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);
1195748696b1SAlp 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);
1196748696b1SAlp 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);
1197748696b1SAlp 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);
1198748696b1SAlp 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);
1199748696b1SAlp 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);
1200748696b1SAlp 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);
1201748696b1SAlp 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);
1202748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1203748696b1SAlp 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);
1204748696b1SAlp 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);
1205748696b1SAlp 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);
1206748696b1SAlp 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);
1207748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1208748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1209748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1210748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1211748696b1SAlp 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);
1212748696b1SAlp 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);
1213c0f10754SAlp 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);
1214eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
121533c78596SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1216eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1217eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1218eb910715SAlp Dener   PetscFunctionReturn(0);
1219eb910715SAlp Dener }
1220eb910715SAlp Dener 
1221eb910715SAlp Dener /*------------------------------------------------------------*/
1222df278d8fSAlp Dener 
1223e0ed867bSAlp Dener PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1224eb910715SAlp Dener {
1225eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1226eb910715SAlp Dener   PetscInt       nrejects;
1227eb910715SAlp Dener   PetscBool      isascii;
1228eb910715SAlp Dener   PetscErrorCode ierr;
1229eb910715SAlp Dener 
1230eb910715SAlp Dener   PetscFunctionBegin;
1231eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1232eb910715SAlp Dener   if (isascii) {
1233eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1234b9ac7092SAlp Dener     if (bnk->M) {
1235cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(bnk->M,&nrejects);CHKERRQ(ierr);
1236b9ac7092SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %D\n",nrejects);CHKERRQ(ierr);
1237eb910715SAlp Dener     }
1238e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1239eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1240b9ac7092SAlp Dener     if (bnk->M) {
1241eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1242b9ac7092SAlp Dener     }
1243eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1244eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1245eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1246eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1247eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1248eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1249eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1250eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1251eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1252eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1253eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1254eb910715SAlp Dener   }
1255eb910715SAlp Dener   PetscFunctionReturn(0);
1256eb910715SAlp Dener }
1257eb910715SAlp Dener 
1258eb910715SAlp Dener /* ---------------------------------------------------------- */
1259df278d8fSAlp Dener 
1260eb910715SAlp Dener /*MC
1261eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
126266ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1263eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1264eb910715SAlp Dener               Hk dk = -gk
12652b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12662b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1267eb910715SAlp Dener 
1268eb910715SAlp Dener     Options Database Keys:
1269e0ed867bSAlp Dener + -max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1270e0ed867bSAlp Dener . -init_type - trust radius initialization method ("constant", "direction", "interpolation")
1271e0ed867bSAlp Dener . -update_type - trust radius update method ("step", "direction", "interpolation")
1272e0ed867bSAlp Dener . -as_type - active-set estimation method ("none", "bertsekas")
1273e0ed867bSAlp Dener . -as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
1274e0ed867bSAlp Dener . -as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
1275e0ed867bSAlp Dener . -sval - (developer) Hessian perturbation starting value
1276e0ed867bSAlp Dener . -imin - (developer) minimum initial Hessian perturbation
1277e0ed867bSAlp Dener . -imax - (developer) maximum initial Hessian perturbation
1278e0ed867bSAlp Dener . -pmin - (developer) minimum Hessian perturbation
1279e0ed867bSAlp Dener . -pmax - (developer) aximum Hessian perturbation
1280e0ed867bSAlp Dener . -pgfac - (developer) Hessian perturbation growth factor
1281e0ed867bSAlp Dener . -psfac - (developer) Hessian perturbation shrink factor
1282e0ed867bSAlp Dener . -imfac - (developer) initial merit factor for Hessian perturbation
1283e0ed867bSAlp Dener . -pmgfac - (developer) merit growth factor for Hessian perturbation
1284e0ed867bSAlp Dener . -pmsfac - (developer) merit shrink factor for Hessian perturbation
1285e0ed867bSAlp Dener . -eta1 - (developer) threshold for rejecting step (-update_type reduction)
1286e0ed867bSAlp Dener . -eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
1287e0ed867bSAlp Dener . -eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
1288e0ed867bSAlp Dener . -eta4 - (developer) threshold for accepting good step (-update_type reduction)
1289e0ed867bSAlp Dener . -alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
1290e0ed867bSAlp Dener . -alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
1291e0ed867bSAlp Dener . -alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
1292e0ed867bSAlp Dener . -alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
1293e0ed867bSAlp Dener . -alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
1294e0ed867bSAlp Dener . -epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
1295e0ed867bSAlp Dener . -mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
1296e0ed867bSAlp Dener . -mu2 - (developer) threshold for accepting good step (-update_type interpolation)
1297e0ed867bSAlp Dener . -gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
1298e0ed867bSAlp Dener . -gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
1299e0ed867bSAlp Dener . -gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
1300e0ed867bSAlp Dener . -gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
1301e0ed867bSAlp Dener . -theta - (developer) trust region interpolation factor (-update_type interpolation)
1302e0ed867bSAlp Dener . -nu1 - (developer) threshold for small line-search step length (-update_type step)
1303e0ed867bSAlp Dener . -nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
1304e0ed867bSAlp Dener . -nu3 - (developer) threshold for large line-search step length (-update_type step)
1305e0ed867bSAlp Dener . -nu4 - (developer) threshold for very large line-search step length (-update_type step)
1306e0ed867bSAlp Dener . -omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
1307e0ed867bSAlp Dener . -omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
1308e0ed867bSAlp Dener . -omega3 - (developer) radius factor for decent line-search step length (-update_type step)
1309e0ed867bSAlp Dener . -omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
1310e0ed867bSAlp Dener . -omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
1311e0ed867bSAlp Dener . -mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
1312e0ed867bSAlp Dener . -mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
1313e0ed867bSAlp Dener . -gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
1314e0ed867bSAlp Dener . -gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
1315e0ed867bSAlp Dener . -gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
1316e0ed867bSAlp Dener . -gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
1317e0ed867bSAlp Dener - -theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1318eb910715SAlp Dener 
1319eb910715SAlp Dener   Level: beginner
1320eb910715SAlp Dener M*/
1321eb910715SAlp Dener 
1322eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1323eb910715SAlp Dener {
1324eb910715SAlp Dener   TAO_BNK        *bnk;
1325eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1326eb910715SAlp Dener   PetscErrorCode ierr;
1327b9ac7092SAlp Dener   PC             pc;
1328eb910715SAlp Dener 
1329eb910715SAlp Dener   PetscFunctionBegin;
1330eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1331eb910715SAlp Dener 
1332eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1333eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1334eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1335eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1336eb910715SAlp Dener 
1337eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1338eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1339eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1340eb910715SAlp Dener 
1341eb910715SAlp Dener   tao->data = (void*)bnk;
1342eb910715SAlp Dener 
134366ed3702SAlp Dener   /*  Hessian shifting parameters */
1344e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1345e0ed867bSAlp Dener   bnk->computestep = TaoBNKComputeStep;
1346e0ed867bSAlp Dener 
1347eb910715SAlp Dener   bnk->sval   = 0.0;
1348eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1349eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1350eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1351eb910715SAlp Dener 
1352eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1353eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1354eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1355eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1356eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1357eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1358eb910715SAlp Dener 
1359eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1360eb910715SAlp Dener   bnk->nu1 = 0.25;
1361eb910715SAlp Dener   bnk->nu2 = 0.50;
1362eb910715SAlp Dener   bnk->nu3 = 1.00;
1363eb910715SAlp Dener   bnk->nu4 = 1.25;
1364eb910715SAlp Dener 
1365eb910715SAlp Dener   bnk->omega1 = 0.25;
1366eb910715SAlp Dener   bnk->omega2 = 0.50;
1367eb910715SAlp Dener   bnk->omega3 = 1.00;
1368eb910715SAlp Dener   bnk->omega4 = 2.00;
1369eb910715SAlp Dener   bnk->omega5 = 4.00;
1370eb910715SAlp Dener 
1371eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1372eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1373eb910715SAlp Dener   bnk->eta2 = 0.25;
1374eb910715SAlp Dener   bnk->eta3 = 0.50;
1375eb910715SAlp Dener   bnk->eta4 = 0.90;
1376eb910715SAlp Dener 
1377eb910715SAlp Dener   bnk->alpha1 = 0.25;
1378eb910715SAlp Dener   bnk->alpha2 = 0.50;
1379eb910715SAlp Dener   bnk->alpha3 = 1.00;
1380eb910715SAlp Dener   bnk->alpha4 = 2.00;
1381eb910715SAlp Dener   bnk->alpha5 = 4.00;
1382eb910715SAlp Dener 
1383eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1384eb910715SAlp Dener   bnk->mu1 = 0.10;
1385eb910715SAlp Dener   bnk->mu2 = 0.50;
1386eb910715SAlp Dener 
1387eb910715SAlp Dener   bnk->gamma1 = 0.25;
1388eb910715SAlp Dener   bnk->gamma2 = 0.50;
1389eb910715SAlp Dener   bnk->gamma3 = 2.00;
1390eb910715SAlp Dener   bnk->gamma4 = 4.00;
1391eb910715SAlp Dener 
1392eb910715SAlp Dener   bnk->theta = 0.05;
1393eb910715SAlp Dener 
1394eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1395eb910715SAlp Dener   bnk->mu1_i = 0.35;
1396eb910715SAlp Dener   bnk->mu2_i = 0.50;
1397eb910715SAlp Dener 
1398eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1399eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1400eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1401eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1402eb910715SAlp Dener 
1403eb910715SAlp Dener   bnk->theta_i = 0.25;
1404eb910715SAlp Dener 
1405eb910715SAlp Dener   /*  Remaining parameters */
1406c0f10754SAlp Dener   bnk->max_cg_its = 0;
1407eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1408eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1409770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14100a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14110a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
141262675beeSAlp Dener   bnk->dmin = 1.0e-6;
141362675beeSAlp Dener   bnk->dmax = 1.0e6;
1414eb910715SAlp Dener 
141583c8fe1dSLisandro Dalcin   bnk->M               = NULL;
141683c8fe1dSLisandro Dalcin   bnk->bfgs_pre        = NULL;
1417eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
14187b1c7716SAlp Dener   bnk->update_type     = BNK_UPDATE_REDUCTION;
14192f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1420eb910715SAlp Dener 
14215eb5f4d6SAlp Dener   bnk->is_stcg         = PETSC_FALSE;
14225eb5f4d6SAlp Dener   bnk->is_gltr         = PETSC_FALSE;
14235eb5f4d6SAlp Dener   bnk->is_nash         = PETSC_FALSE;
14245eb5f4d6SAlp Dener 
1425e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1426e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1427e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1428e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1429e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1430e031d6f5SAlp Dener 
1431c0f10754SAlp Dener   /* Create the line search */
1432eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1433eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1434e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1435eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1436eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1437eb910715SAlp Dener 
1438eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1439eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1440eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1441e0ed867bSAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,"tao_bnk_");CHKERRQ(ierr);
144205de396fSBarry Smith   ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr);
1443f5a7d1c1SBarry Smith   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
1444b9ac7092SAlp Dener   ierr = PCSetType(pc, PCLMVM);CHKERRQ(ierr);
1445eb910715SAlp Dener   PetscFunctionReturn(0);
1446eb910715SAlp Dener }
1447