xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 4f4fdda464867e01e2f7b514a4f52c421e99e383)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener #include <petscksp.h>
5eb910715SAlp Dener 
6e031d6f5SAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
7e031d6f5SAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
8e031d6f5SAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
9e031d6f5SAlp Dener 
10e031d6f5SAlp Dener /*------------------------------------------------------------*/
11e031d6f5SAlp Dener 
12df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
13df278d8fSAlp Dener 
14c0f10754SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
15eb910715SAlp Dener {
16eb910715SAlp Dener   PetscErrorCode               ierr;
17eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
18eb910715SAlp Dener   PC                           pc;
19eb910715SAlp Dener 
2089da521bSAlp Dener   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma, resnorm;
21eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
220ad3a497SAlp Dener   PetscBool                    is_bfgs, is_jacobi, is_symmetric, sym_set;
23c4b75bccSAlp Dener   PetscInt                     n, N, nDiff;
24eb910715SAlp Dener   PetscInt                     i_max = 5;
25eb910715SAlp Dener   PetscInt                     j_max = 1;
26eb910715SAlp Dener   PetscInt                     i, j;
27eb910715SAlp Dener 
28eb910715SAlp Dener   PetscFunctionBegin;
2928017e9fSAlp Dener   /* Project the current point onto the feasible set */
3028017e9fSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
31e031d6f5SAlp Dener   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
32b9ac7092SAlp Dener   if (tao->bounded) {
3328017e9fSAlp Dener     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
34cd929ea3SAlp Dener   }
3528017e9fSAlp Dener 
3628017e9fSAlp Dener   /* Project the initial point onto the feasible region */
373b063059SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
3828017e9fSAlp Dener 
3928017e9fSAlp Dener   /* Check convergence criteria */
4028017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
4161be54a6SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
4261be54a6SAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
4361be54a6SAlp Dener   ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
4408752603SAlp Dener   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
4528017e9fSAlp Dener 
46c0f10754SAlp Dener   /* Test the initial point for convergence */
4789da521bSAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
4889da521bSAlp Dener   ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
49b4a30f08SAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
50e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
51e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
52c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
53c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
54c0f10754SAlp Dener 
55e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
56eb910715SAlp Dener   bnk->ksp_atol = 0;
57eb910715SAlp Dener   bnk->ksp_rtol = 0;
58eb910715SAlp Dener   bnk->ksp_dtol = 0;
59eb910715SAlp Dener   bnk->ksp_ctol = 0;
60eb910715SAlp Dener   bnk->ksp_negc = 0;
61eb910715SAlp Dener   bnk->ksp_iter = 0;
62eb910715SAlp Dener   bnk->ksp_othr = 0;
63eb910715SAlp Dener 
64e031d6f5SAlp Dener   /* Reset accepted step type counters */
65e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
66e031d6f5SAlp Dener   bnk->newt = 0;
67e031d6f5SAlp Dener   bnk->bfgs = 0;
68e031d6f5SAlp Dener   bnk->sgrad = 0;
69e031d6f5SAlp Dener   bnk->grad = 0;
70e031d6f5SAlp Dener 
71fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
72fed79b8eSAlp Dener   bnk->pert = bnk->sval;
73fed79b8eSAlp Dener 
74937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
75937a31a1SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
76937a31a1SAlp Dener 
77e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
78b9ac7092SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
79b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr);
80b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr);
81b9ac7092SAlp Dener   if (is_bfgs) {
82b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
83b9ac7092SAlp Dener     ierr = PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M);CHKERRQ(ierr);
84eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
85eb910715SAlp Dener     ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
86b9ac7092SAlp Dener     ierr = MatSetSizes(bnk->M, n, n, N, N);CHKERRQ(ierr);
87cd929ea3SAlp Dener     ierr = MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
880ad3a497SAlp Dener     ierr = MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric);CHKERRQ(ierr);
890ad3a497SAlp Dener     if (!sym_set || !is_symmetric) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
90b9ac7092SAlp Dener   } else if (is_jacobi) {
91b9ac7092SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
92e031d6f5SAlp Dener   }
93e031d6f5SAlp Dener 
94e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
9562675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
9662675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
97eb910715SAlp Dener 
98eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
99eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
100c0f10754SAlp Dener   *needH = PETSC_TRUE;
101eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
10262675beeSAlp Dener     switch(initType) {
103eb910715SAlp Dener     case BNK_INIT_CONSTANT:
104eb910715SAlp Dener       /* Use the initial radius specified */
105c0f10754SAlp Dener       tao->trust = tao->trust0;
106eb910715SAlp Dener       break;
107eb910715SAlp Dener 
108eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
109c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
110eb910715SAlp Dener       max_radius = 0.0;
11108752603SAlp Dener       tao->trust = tao->trust0;
112eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1130a4511e9SAlp Dener         f_min = bnk->f;
114eb910715SAlp Dener         sigma = 0.0;
115eb910715SAlp Dener 
116c0f10754SAlp Dener         if (*needH) {
11762602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
118e0ed867bSAlp Dener           ierr = bnk->computehessian(tao);CHKERRQ(ierr);
11908752603SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
12089da521bSAlp Dener           ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
12189da521bSAlp Dener           if (bnk->active_idx) {
1222ab2a32cSAlp Dener             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
12328017e9fSAlp Dener           } else {
12408752603SAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
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);
139b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "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 
176eb910715SAlp Dener             if (PetscAbsScalar(kappa - 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               }
1878f8a4e06SAlp Dener             } else if (PetscAbsScalar(kappa - 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 */
231c4b75bccSAlp Dener           ierr = VecNorm(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);
236b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "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 */
242937a31a1SAlp Dener           ierr = TaoBNCGSetRecycleFlag(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 */
278e0ed867bSAlp Dener   if (bnk->active_idx) {
279e0ed867bSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
280e0ed867bSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
281e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
282e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
283e0ed867bSAlp Dener     } else {
284e0ed867bSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
285e0ed867bSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
286e0ed867bSAlp Dener     }
287e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
288e0ed867bSAlp Dener       ierr = PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx);CHKERRQ(ierr);
289e0ed867bSAlp Dener     }
290e0ed867bSAlp Dener   } else {
291e0ed867bSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
292e0ed867bSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
293e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
294e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
295e0ed867bSAlp Dener     } else {
296e0ed867bSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
297e0ed867bSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);CHKERRQ(ierr);
298e0ed867bSAlp Dener     }
299e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
300e0ed867bSAlp Dener       ierr = PCLMVMClearIS(bnk->bfgs_pre);CHKERRQ(ierr);
301e0ed867bSAlp Dener     }
302e0ed867bSAlp Dener   }
30362675beeSAlp Dener   PetscFunctionReturn(0);
30462675beeSAlp Dener }
30562675beeSAlp Dener 
30662675beeSAlp Dener /*------------------------------------------------------------*/
30762675beeSAlp Dener 
3082f75a4aaSAlp Dener /* Routine for estimating the active set */
3092f75a4aaSAlp Dener 
31008752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3112f75a4aaSAlp Dener {
3122f75a4aaSAlp Dener   PetscErrorCode               ierr;
3132f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
314c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3152f75a4aaSAlp Dener 
3162f75a4aaSAlp Dener   PetscFunctionBegin;
31708752603SAlp Dener   switch (asType) {
3182f75a4aaSAlp Dener   case BNK_AS_NONE:
3192f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3202f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3212f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3222f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3232f75a4aaSAlp Dener     break;
3242f75a4aaSAlp Dener 
3252f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3262f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
327b9ac7092SAlp Dener     if (bnk->M) {
3282f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
329cd929ea3SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3302f75a4aaSAlp Dener     } else {
33161be54a6SAlp Dener       ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
332c4b75bccSAlp Dener       ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
333c4b75bccSAlp Dener       if (hessComputed && diagExists) {
3349b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3352f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3369b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3379b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3382f75a4aaSAlp Dener         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3392f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
34061be54a6SAlp Dener       } else {
341c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
34261be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
34361be54a6SAlp Dener       }
3442f75a4aaSAlp Dener     }
3452f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
34689da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
34789da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
348c4b75bccSAlp Dener     break;
3492f75a4aaSAlp Dener 
3502f75a4aaSAlp Dener   default:
3512f75a4aaSAlp Dener     break;
3522f75a4aaSAlp Dener   }
3532f75a4aaSAlp Dener   PetscFunctionReturn(0);
3542f75a4aaSAlp Dener }
3552f75a4aaSAlp Dener 
3562f75a4aaSAlp Dener /*------------------------------------------------------------*/
3572f75a4aaSAlp Dener 
3582f75a4aaSAlp Dener /* Routine for bounding the step direction */
3592f75a4aaSAlp Dener 
360a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
3612f75a4aaSAlp Dener {
3622f75a4aaSAlp Dener   PetscErrorCode               ierr;
3632f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3642f75a4aaSAlp Dener 
3652f75a4aaSAlp Dener   PetscFunctionBegin;
366a1318120SAlp Dener   switch (asType) {
3672f75a4aaSAlp Dener   case BNK_AS_NONE:
368c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
3692f75a4aaSAlp Dener     break;
3702f75a4aaSAlp Dener 
3712f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
372c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
3732f75a4aaSAlp Dener     break;
3742f75a4aaSAlp Dener 
3752f75a4aaSAlp Dener   default:
3762f75a4aaSAlp Dener     break;
3772f75a4aaSAlp Dener   }
3782f75a4aaSAlp Dener   PetscFunctionReturn(0);
3792f75a4aaSAlp Dener }
3802f75a4aaSAlp Dener 
381e031d6f5SAlp Dener /*------------------------------------------------------------*/
382e031d6f5SAlp Dener 
383e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
384e031d6f5SAlp Dener    accelerate Newton convergence.
385e031d6f5SAlp Dener 
386e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
387e031d6f5SAlp Dener    for more gradient evaluations.
388e031d6f5SAlp Dener */
389e031d6f5SAlp Dener 
390c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
391c0f10754SAlp Dener {
392c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
393c0f10754SAlp Dener   PetscErrorCode               ierr;
394c0f10754SAlp Dener 
395c0f10754SAlp Dener   PetscFunctionBegin;
396c0f10754SAlp Dener   *terminate = PETSC_FALSE;
397c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
398c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
399c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
400c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
401c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
402c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
403c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
404c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
405c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
406c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
407e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
408c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
409c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
410c0f10754SAlp 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) {
411c0f10754SAlp Dener       *terminate = PETSC_TRUE;
41261be54a6SAlp Dener     } else {
41333c78596SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
414c0f10754SAlp Dener     }
415c0f10754SAlp Dener   }
416c0f10754SAlp Dener   PetscFunctionReturn(0);
417c0f10754SAlp Dener }
418c0f10754SAlp Dener 
4192f75a4aaSAlp Dener /*------------------------------------------------------------*/
4202f75a4aaSAlp Dener 
421c0f10754SAlp Dener /* Routine for computing the Newton step. */
422df278d8fSAlp Dener 
42362675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
424eb910715SAlp Dener {
425eb910715SAlp Dener   PetscErrorCode               ierr;
426eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
427eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
428eb910715SAlp Dener   PetscInt                     kspits;
429eb910715SAlp Dener 
430eb910715SAlp Dener   PetscFunctionBegin;
43189da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
43289da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
43389da521bSAlp Dener   if (!bnk->inactive_idx) {
43489da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
435a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
43689da521bSAlp Dener     PetscFunctionReturn(0);
43789da521bSAlp Dener   }
43889da521bSAlp Dener 
43962675beeSAlp Dener   /* Shift the reduced Hessian matrix */
44062675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
44162675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
44262675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
44362675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
44462675beeSAlp Dener     }
44562675beeSAlp Dener   }
44662675beeSAlp Dener 
447eb910715SAlp Dener   /* Solve the Newton system of equations */
448937a31a1SAlp Dener   tao->ksp_its = 0;
4492f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
4505e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
45109164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
4525e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
45389da521bSAlp Dener   if (bnk->active_idx) {
4545e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
4555e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
4565e9b73cbSAlp Dener   } else {
4575e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4585e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
45928017e9fSAlp Dener   }
460eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
461fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4625e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
463eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
464eb910715SAlp Dener     tao->ksp_its+=kspits;
465eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
466080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
467eb910715SAlp Dener 
468eb910715SAlp Dener     if (0.0 == tao->trust) {
469eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
470080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
471080d2917SAlp Dener         tao->trust = bnk->dnorm;
472eb910715SAlp Dener 
473eb910715SAlp Dener         /* Modify the radius if it is too large or small */
474eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
475eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
476eb910715SAlp Dener       } else {
477eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
478eb910715SAlp Dener            the trust-region subproblem to get a direction */
479eb910715SAlp Dener         tao->trust = tao->trust0;
480eb910715SAlp Dener 
481eb910715SAlp Dener         /* Modify the radius if it is too large or small */
482eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
483eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
484eb910715SAlp Dener 
485fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4865e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
487eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
488eb910715SAlp Dener         tao->ksp_its+=kspits;
489eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
490080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
491eb910715SAlp Dener 
492080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
493eb910715SAlp Dener       }
494eb910715SAlp Dener     }
495eb910715SAlp Dener   } else {
4965e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
497eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
498eb910715SAlp Dener     tao->ksp_its += kspits;
499eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
500eb910715SAlp Dener   }
5015e9b73cbSAlp Dener   /* Restore sub vectors back */
50289da521bSAlp Dener   if (bnk->active_idx) {
5035e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5045e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5055e9b73cbSAlp Dener   }
506770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
507fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
508a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
509770b7498SAlp Dener 
510770b7498SAlp Dener   /* Record convergence reasons */
511e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
512e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
513770b7498SAlp Dener     ++bnk->ksp_atol;
514e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
515770b7498SAlp Dener     ++bnk->ksp_rtol;
516e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
517770b7498SAlp Dener     ++bnk->ksp_ctol;
518e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
519770b7498SAlp Dener     ++bnk->ksp_negc;
520e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
521770b7498SAlp Dener     ++bnk->ksp_dtol;
522e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
523770b7498SAlp Dener     ++bnk->ksp_iter;
524770b7498SAlp Dener   } else {
525770b7498SAlp Dener     ++bnk->ksp_othr;
526770b7498SAlp Dener   }
527fed79b8eSAlp Dener 
528fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
529b9ac7092SAlp Dener   if (bnk->M) {
530cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
531b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
532fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
533cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
53409164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
535eb910715SAlp Dener     }
536fed79b8eSAlp Dener   }
537e465cd6fSAlp Dener   PetscFunctionReturn(0);
538e465cd6fSAlp Dener }
539eb910715SAlp Dener 
54062675beeSAlp Dener /*------------------------------------------------------------*/
54162675beeSAlp Dener 
5425e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5435e9b73cbSAlp Dener 
5445e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
5455e9b73cbSAlp Dener {
5465e9b73cbSAlp Dener   PetscErrorCode               ierr;
5475e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
5485e9b73cbSAlp Dener 
5495e9b73cbSAlp Dener   PetscFunctionBegin;
5505e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
55189da521bSAlp Dener   if (bnk->active_idx){
5525e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5535e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5545e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5555e9b73cbSAlp Dener   } else {
5565e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
5575e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5585e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
5595e9b73cbSAlp Dener   }
5605e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5615e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
5625e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
56333c78596SAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);CHKERRQ(ierr);
5645e9b73cbSAlp Dener   /* Restore the sub vectors */
56589da521bSAlp Dener   if (bnk->active_idx){
5665e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5675e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5685e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5695e9b73cbSAlp Dener   }
5705e9b73cbSAlp Dener   PetscFunctionReturn(0);
5715e9b73cbSAlp Dener }
5725e9b73cbSAlp Dener 
5735e9b73cbSAlp Dener /*------------------------------------------------------------*/
5745e9b73cbSAlp Dener 
57562675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
57662675beeSAlp Dener 
57762675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
57862675beeSAlp Dener    in the event that the Newton step fails the test.
57962675beeSAlp Dener */
58062675beeSAlp Dener 
581e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
582e465cd6fSAlp Dener {
583e465cd6fSAlp Dener   PetscErrorCode               ierr;
584e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
585e465cd6fSAlp Dener 
586b2d8c577SAlp Dener   PetscReal                    gdx, e_min;
587e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
588e465cd6fSAlp Dener 
589e465cd6fSAlp Dener   PetscFunctionBegin;
590080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
591eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
592eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
593eb910715SAlp Dener        Update the perturbation for next time */
594eb910715SAlp Dener     if (bnk->pert <= 0.0) {
595eb910715SAlp Dener       /* Initialize the perturbation */
596eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
597eb910715SAlp Dener       if (bnk->is_gltr) {
598eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
599eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
600eb910715SAlp Dener       }
601eb910715SAlp Dener     } else {
602eb910715SAlp Dener       /* Increase the perturbation */
603eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
604eb910715SAlp Dener     }
605eb910715SAlp Dener 
6060ad3a497SAlp Dener     if (!bnk->M) {
607eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
608eb910715SAlp Dener          Must use gradient direction in this case */
609080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
610eb910715SAlp Dener       *stepType = BNK_GRADIENT;
611eb910715SAlp Dener     } else {
612eb910715SAlp Dener       /* Attempt to use the BFGS direction */
613cd929ea3SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
614eb910715SAlp Dener 
6158d5ead36SAlp Dener       /* Check for success (descent direction)
6168d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
6178d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
618080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
6193105154fSTodd Munson       if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
620eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
621eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
622eb910715SAlp Dener            the first solve produces the scaled gradient direction,
623eb910715SAlp Dener            which is guaranteed to be descent */
624eb910715SAlp Dener 
625eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
626cd929ea3SAlp Dener         ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
62709164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
628cd929ea3SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
629eb910715SAlp Dener 
630eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
631eb910715SAlp Dener       } else {
632cd929ea3SAlp Dener         ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
633eb910715SAlp Dener         if (1 == bfgsUpdates) {
634eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
635eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
636eb910715SAlp Dener         } else {
637eb910715SAlp Dener           *stepType = BNK_BFGS;
638eb910715SAlp Dener         }
639eb910715SAlp Dener       }
640eb910715SAlp Dener     }
6418d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6428d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
643a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
644eb910715SAlp Dener   } else {
645eb910715SAlp Dener     /* Computed Newton step is descent */
646eb910715SAlp Dener     switch (ksp_reason) {
647eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
648eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
649eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
650eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
651eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
652eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
653eb910715SAlp Dener       if (bnk->pert <= 0.0) {
654eb910715SAlp Dener         /* Initialize the perturbation */
655eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
656eb910715SAlp Dener         if (bnk->is_gltr) {
657eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
658eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
659eb910715SAlp Dener         }
660eb910715SAlp Dener       } else {
661eb910715SAlp Dener         /* Increase the perturbation */
662eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
663eb910715SAlp Dener       }
664eb910715SAlp Dener       break;
665eb910715SAlp Dener 
666eb910715SAlp Dener     default:
667eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
668eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
669eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
670eb910715SAlp Dener         bnk->pert = 0.0;
671eb910715SAlp Dener       }
672eb910715SAlp Dener       break;
673eb910715SAlp Dener     }
674fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
675eb910715SAlp Dener   }
676eb910715SAlp Dener   PetscFunctionReturn(0);
677eb910715SAlp Dener }
678eb910715SAlp Dener 
679df278d8fSAlp Dener /*------------------------------------------------------------*/
680df278d8fSAlp Dener 
681df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
682df278d8fSAlp Dener 
683df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
684df278d8fSAlp Dener   Newton step does not produce a valid step length.
685df278d8fSAlp Dener */
686df278d8fSAlp Dener 
687937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
688c14b763aSAlp Dener {
689c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
690c14b763aSAlp Dener   PetscErrorCode ierr;
691c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
692c14b763aSAlp Dener 
693b2d8c577SAlp Dener   PetscReal      e_min, gdx;
694c14b763aSAlp Dener   PetscInt       bfgsUpdates;
695c14b763aSAlp Dener 
696c14b763aSAlp Dener   PetscFunctionBegin;
697c14b763aSAlp Dener   /* Perform the linesearch */
698c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
699c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
700c14b763aSAlp Dener 
701b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
702c14b763aSAlp Dener     /* Linesearch failed, revert solution */
703c14b763aSAlp Dener     bnk->f = bnk->fold;
704c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
705c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
706c14b763aSAlp Dener 
707937a31a1SAlp Dener     switch(*stepType) {
708c14b763aSAlp Dener     case BNK_NEWTON:
7098d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
710c14b763aSAlp Dener          Update the perturbation for next time */
711c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
712c14b763aSAlp Dener         /* Initialize the perturbation */
713c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
714c14b763aSAlp Dener         if (bnk->is_gltr) {
715c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
716c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
717c14b763aSAlp Dener         }
718c14b763aSAlp Dener       } else {
719c14b763aSAlp Dener         /* Increase the perturbation */
720c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
721c14b763aSAlp Dener       }
722c14b763aSAlp Dener 
7230ad3a497SAlp Dener       if (!bnk->M) {
724c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
725c14b763aSAlp Dener            Must use gradient direction in this case */
726937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
727937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
728c14b763aSAlp Dener       } else {
729c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
730cd929ea3SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
7318d5ead36SAlp Dener         /* Check for success (descent direction)
7328d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
733c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7343105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
735c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
736c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
737c14b763aSAlp Dener              Use steepest descent direction (scaled) */
738cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
739c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
740cd929ea3SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
741c14b763aSAlp Dener 
742c14b763aSAlp Dener           bfgsUpdates = 1;
743937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
744c14b763aSAlp Dener         } else {
745cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
746c14b763aSAlp Dener           if (1 == bfgsUpdates) {
747c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
748937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
749c14b763aSAlp Dener           } else {
750937a31a1SAlp Dener             *stepType = BNK_BFGS;
751c14b763aSAlp Dener           }
752c14b763aSAlp Dener         }
753c14b763aSAlp Dener       }
754c14b763aSAlp Dener       break;
755c14b763aSAlp Dener 
756c14b763aSAlp Dener     case BNK_BFGS:
757c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
758c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
759c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
760cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
761c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
762cd929ea3SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
763c14b763aSAlp Dener 
764c14b763aSAlp Dener       bfgsUpdates = 1;
765937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
766c14b763aSAlp Dener       break;
767c14b763aSAlp Dener     }
7688d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7698d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
770a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
771c14b763aSAlp Dener 
7728d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
773c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
774c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
775c14b763aSAlp Dener   }
776c14b763aSAlp Dener   *reason = ls_reason;
777c14b763aSAlp Dener   PetscFunctionReturn(0);
778c14b763aSAlp Dener }
779c14b763aSAlp Dener 
780df278d8fSAlp Dener /*------------------------------------------------------------*/
781df278d8fSAlp Dener 
782df278d8fSAlp Dener /* Routine for updating the trust radius.
783df278d8fSAlp Dener 
784df278d8fSAlp Dener   Function features three different update methods:
785df278d8fSAlp Dener   1) Line-search step length based
786df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
787df278d8fSAlp Dener   3) Interpolation
788df278d8fSAlp Dener */
789df278d8fSAlp Dener 
79028017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
791080d2917SAlp Dener {
792080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
793080d2917SAlp Dener   PetscErrorCode ierr;
794080d2917SAlp Dener 
795b1c2d0e3SAlp Dener   PetscReal      step, kappa;
796080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
797080d2917SAlp Dener 
798080d2917SAlp Dener   PetscFunctionBegin;
799080d2917SAlp Dener   /* Update trust region radius */
800080d2917SAlp Dener   *accept = PETSC_FALSE;
80128017e9fSAlp Dener   switch(updateType) {
802080d2917SAlp Dener   case BNK_UPDATE_STEP:
803c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
804080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
805080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
806080d2917SAlp Dener       if (step < bnk->nu1) {
807080d2917SAlp Dener         /* Very bad step taken; reduce radius */
808080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
809080d2917SAlp Dener       } else if (step < bnk->nu2) {
810080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
811080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
812080d2917SAlp Dener       } else if (step < bnk->nu3) {
813080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
814080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
815080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
816080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
817080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
818080d2917SAlp Dener         }
819080d2917SAlp Dener       } else if (step < bnk->nu4) {
820080d2917SAlp Dener         /*  Full step taken; increase the radius */
821080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
822080d2917SAlp Dener       } else {
823080d2917SAlp Dener         /*  More than full step taken; increase the radius */
824080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
825080d2917SAlp Dener       }
826080d2917SAlp Dener     } else {
827080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
828080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
829080d2917SAlp Dener     }
830080d2917SAlp Dener     break;
831080d2917SAlp Dener 
832080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
833080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
834e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
835fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
836fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
837fed79b8eSAlp Dener            be rejected! */
838080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8393105154fSTodd Munson       } else {
840b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
841080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
842080d2917SAlp Dener         } else {
8433105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
844080d2917SAlp Dener             kappa = 1.0;
8453105154fSTodd Munson           } else {
846080d2917SAlp Dener             kappa = actred / prered;
847080d2917SAlp Dener           }
848fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
849080d2917SAlp Dener           if (kappa < bnk->eta1) {
850fed79b8eSAlp Dener             /* Reject the step */
851080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8523105154fSTodd Munson           } else {
853fed79b8eSAlp Dener             /* Accept the step */
854c133c014SAlp Dener             *accept = PETSC_TRUE;
855c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8568d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
857080d2917SAlp Dener               if (kappa < bnk->eta2) {
858080d2917SAlp Dener                 /* Marginal bad step */
859c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
8603105154fSTodd Munson               } else if (kappa < bnk->eta3) {
861fed79b8eSAlp Dener                 /* Reasonable step */
862fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
8633105154fSTodd Munson               } else if (kappa < bnk->eta4) {
864080d2917SAlp Dener                 /* Good step */
865c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
8663105154fSTodd Munson               } else {
867080d2917SAlp Dener                 /* Very good step */
868c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
869080d2917SAlp Dener               }
870c133c014SAlp Dener             }
871080d2917SAlp Dener           }
872080d2917SAlp Dener         }
873080d2917SAlp Dener       }
874080d2917SAlp Dener     } else {
875080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
876080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
877080d2917SAlp Dener     }
878080d2917SAlp Dener     break;
879080d2917SAlp Dener 
880080d2917SAlp Dener   default:
881080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
882b1c2d0e3SAlp Dener       if (prered < 0.0) {
883080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
884080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
885080d2917SAlp Dener         /*  be rejected! */
886080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
887080d2917SAlp Dener       } else {
888b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
889080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
890080d2917SAlp Dener         } else {
891080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
892080d2917SAlp Dener             kappa = 1.0;
893080d2917SAlp Dener           } else {
894080d2917SAlp Dener             kappa = actred / prered;
895080d2917SAlp Dener           }
896080d2917SAlp Dener 
897080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
898080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
899080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
900080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
901080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
902080d2917SAlp Dener 
903080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
904080d2917SAlp Dener             /*  Great agreement */
905080d2917SAlp Dener             *accept = PETSC_TRUE;
906080d2917SAlp Dener             if (tau_max < 1.0) {
907080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
908080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
909080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
910080d2917SAlp Dener             } else {
911080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
912080d2917SAlp Dener             }
913080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
914080d2917SAlp Dener             /*  Good agreement */
915080d2917SAlp Dener             *accept = PETSC_TRUE;
916080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
917080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
918080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
919080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
920080d2917SAlp Dener             } else if (tau_max < 1.0) {
921080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
922080d2917SAlp Dener             } else {
923080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
924080d2917SAlp Dener             }
925080d2917SAlp Dener           } else {
926080d2917SAlp Dener             /*  Not good agreement */
927080d2917SAlp Dener             if (tau_min > 1.0) {
928080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
929080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
930080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
931080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
932080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
933080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
934080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
935080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
936080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
937080d2917SAlp Dener             } else {
938080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
939080d2917SAlp Dener             }
940080d2917SAlp Dener           }
941080d2917SAlp Dener         }
942080d2917SAlp Dener       }
943080d2917SAlp Dener     } else {
944080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
945080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
946080d2917SAlp Dener     }
94728017e9fSAlp Dener     break;
948080d2917SAlp Dener   }
949c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
950c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
951fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
952080d2917SAlp Dener   PetscFunctionReturn(0);
953080d2917SAlp Dener }
954080d2917SAlp Dener 
955eb910715SAlp Dener /* ---------------------------------------------------------- */
956df278d8fSAlp Dener 
95762675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
95862675beeSAlp Dener {
95962675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
96062675beeSAlp Dener 
96162675beeSAlp Dener   PetscFunctionBegin;
96262675beeSAlp Dener   switch (stepType) {
96362675beeSAlp Dener   case BNK_NEWTON:
96462675beeSAlp Dener     ++bnk->newt;
96562675beeSAlp Dener     break;
96662675beeSAlp Dener   case BNK_BFGS:
96762675beeSAlp Dener     ++bnk->bfgs;
96862675beeSAlp Dener     break;
96962675beeSAlp Dener   case BNK_SCALED_GRADIENT:
97062675beeSAlp Dener     ++bnk->sgrad;
97162675beeSAlp Dener     break;
97262675beeSAlp Dener   case BNK_GRADIENT:
97362675beeSAlp Dener     ++bnk->grad;
97462675beeSAlp Dener     break;
97562675beeSAlp Dener   default:
97662675beeSAlp Dener     break;
97762675beeSAlp Dener   }
97862675beeSAlp Dener   PetscFunctionReturn(0);
97962675beeSAlp Dener }
98062675beeSAlp Dener 
98162675beeSAlp Dener /* ---------------------------------------------------------- */
98262675beeSAlp Dener 
9839b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
984eb910715SAlp Dener {
985eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
986eb910715SAlp Dener   PetscErrorCode ierr;
987e031d6f5SAlp Dener   PetscInt       i;
988eb910715SAlp Dener 
989eb910715SAlp Dener   PetscFunctionBegin;
990c4b75bccSAlp Dener   if (!tao->gradient) {
991c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
992c4b75bccSAlp Dener   }
993c4b75bccSAlp Dener   if (!tao->stepdirection) {
994c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
995c4b75bccSAlp Dener   }
996c4b75bccSAlp Dener   if (!bnk->W) {
997c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
998c4b75bccSAlp Dener   }
999c4b75bccSAlp Dener   if (!bnk->Xold) {
1000c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1001c4b75bccSAlp Dener   }
1002c4b75bccSAlp Dener   if (!bnk->Gold) {
1003c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1004c4b75bccSAlp Dener   }
1005c4b75bccSAlp Dener   if (!bnk->Xwork) {
1006c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1007c4b75bccSAlp Dener   }
1008c4b75bccSAlp Dener   if (!bnk->Gwork) {
1009c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1010c4b75bccSAlp Dener   }
1011c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1012c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1013c4b75bccSAlp Dener   }
1014c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1015c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1016c4b75bccSAlp Dener   }
1017c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1018c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1019c4b75bccSAlp Dener   }
1020c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1021c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1022c4b75bccSAlp Dener   }
1023e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1024c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1025c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
102689da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
102789da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
102889da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1029c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1030c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1031c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1032c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1033c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1034c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1035c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1036c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1037c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1038c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1039c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1040c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1041c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1042c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1043e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1044e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1045e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1046937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1047e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1048e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1049e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1050e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1051c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1052e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1053e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1054e031d6f5SAlp Dener     }
1055e031d6f5SAlp Dener   }
1056c0f10754SAlp Dener   bnk->X_inactive = 0;
1057c0f10754SAlp Dener   bnk->G_inactive = 0;
105862675beeSAlp Dener   bnk->inactive_work = 0;
105962675beeSAlp Dener   bnk->active_work = 0;
106062675beeSAlp Dener   bnk->inactive_idx = 0;
106162675beeSAlp Dener   bnk->active_idx = 0;
106262675beeSAlp Dener   bnk->active_lower = 0;
106362675beeSAlp Dener   bnk->active_upper = 0;
106462675beeSAlp Dener   bnk->active_fixed = 0;
1065eb910715SAlp Dener   bnk->M = 0;
106609164190SAlp Dener   bnk->H_inactive = 0;
106709164190SAlp Dener   bnk->Hpre_inactive = 0;
1068eb910715SAlp Dener   PetscFunctionReturn(0);
1069eb910715SAlp Dener }
1070eb910715SAlp Dener 
1071eb910715SAlp Dener /*------------------------------------------------------------*/
1072df278d8fSAlp Dener 
1073e0ed867bSAlp Dener PetscErrorCode TaoDestroy_BNK(Tao tao)
1074eb910715SAlp Dener {
1075eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1076eb910715SAlp Dener   PetscErrorCode ierr;
1077eb910715SAlp Dener 
1078eb910715SAlp Dener   PetscFunctionBegin;
1079eb910715SAlp Dener   if (tao->setupcalled) {
1080eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1081eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1082eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
108309164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
108409164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
108509164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
108609164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
108762675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
108862675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1089c4b75bccSAlp Dener   }
1090ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1091ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1092ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1093ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1094ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
1095c4b75bccSAlp Dener   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {
1096c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1097c4b75bccSAlp Dener   }
1098c4b75bccSAlp Dener   if (bnk->H_inactive != tao->hessian) {
1099c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1100c4b75bccSAlp Dener   }
1101ca964c49SAlp Dener   ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1102eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1103eb910715SAlp Dener   PetscFunctionReturn(0);
1104eb910715SAlp Dener }
1105eb910715SAlp Dener 
1106eb910715SAlp Dener /*------------------------------------------------------------*/
1107df278d8fSAlp Dener 
1108e0ed867bSAlp Dener PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1109eb910715SAlp Dener {
1110eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1111eb910715SAlp Dener   PetscErrorCode ierr;
1112e0ed867bSAlp Dener   KSPType        ksp_type;
1113eb910715SAlp Dener 
1114eb910715SAlp Dener   PetscFunctionBegin;
1115*4f4fdda4SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
11168d5ead36SAlp Dener   ierr = PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, 0);CHKERRQ(ierr);
11178d5ead36SAlp Dener   ierr = PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);CHKERRQ(ierr);
11182f75a4aaSAlp Dener   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr);
1119748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1120748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1121748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1122748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1123748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1124748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1125748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1126748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1127748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1128748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1129748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1130748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1131748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1132748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1133748696b1SAlp 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);
1134748696b1SAlp 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);
1135748696b1SAlp 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);
1136748696b1SAlp 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);
1137748696b1SAlp 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);
1138748696b1SAlp 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);
1139748696b1SAlp 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);
1140748696b1SAlp 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);
1141748696b1SAlp 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);
1142748696b1SAlp 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);
1143748696b1SAlp 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);
1144748696b1SAlp 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);
1145748696b1SAlp 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);
1146748696b1SAlp 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);
1147748696b1SAlp 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);
1148748696b1SAlp 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);
1149748696b1SAlp 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);
1150748696b1SAlp 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);
1151748696b1SAlp 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);
1152748696b1SAlp 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);
1153748696b1SAlp 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);
1154748696b1SAlp 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);
1155748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1156748696b1SAlp 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);
1157748696b1SAlp 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);
1158748696b1SAlp 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);
1159748696b1SAlp 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);
1160748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1161748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1162748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1163748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1164748696b1SAlp 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);
1165748696b1SAlp 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);
1166c0f10754SAlp 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);
1167eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
116833c78596SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1169eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1170eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1171e0ed867bSAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
1172e0ed867bSAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
1173e0ed867bSAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
1174e0ed867bSAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1175eb910715SAlp Dener   PetscFunctionReturn(0);
1176eb910715SAlp Dener }
1177eb910715SAlp Dener 
1178eb910715SAlp Dener /*------------------------------------------------------------*/
1179df278d8fSAlp Dener 
1180e0ed867bSAlp Dener PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1181eb910715SAlp Dener {
1182eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1183eb910715SAlp Dener   PetscInt       nrejects;
1184eb910715SAlp Dener   PetscBool      isascii;
1185eb910715SAlp Dener   PetscErrorCode ierr;
1186eb910715SAlp Dener 
1187eb910715SAlp Dener   PetscFunctionBegin;
1188eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1189eb910715SAlp Dener   if (isascii) {
1190eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1191b9ac7092SAlp Dener     if (bnk->M) {
1192cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(bnk->M,&nrejects);CHKERRQ(ierr);
1193b9ac7092SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %D\n",nrejects);CHKERRQ(ierr);
1194eb910715SAlp Dener     }
1195e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1196eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1197b9ac7092SAlp Dener     if (bnk->M) {
1198eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1199b9ac7092SAlp Dener     }
1200eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1201eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1202eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1203eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1204eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1205eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1206eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1207eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1208eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1209eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1210eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1211eb910715SAlp Dener   }
1212eb910715SAlp Dener   PetscFunctionReturn(0);
1213eb910715SAlp Dener }
1214eb910715SAlp Dener 
1215eb910715SAlp Dener /* ---------------------------------------------------------- */
1216df278d8fSAlp Dener 
1217eb910715SAlp Dener /*MC
1218eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
121966ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1220eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1221eb910715SAlp Dener               Hk dk = -gk
12222b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12232b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1224eb910715SAlp Dener 
1225eb910715SAlp Dener     Options Database Keys:
1226e0ed867bSAlp Dener + -max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1227e0ed867bSAlp Dener . -init_type - trust radius initialization method ("constant", "direction", "interpolation")
1228e0ed867bSAlp Dener . -update_type - trust radius update method ("step", "direction", "interpolation")
1229e0ed867bSAlp Dener . -as_type - active-set estimation method ("none", "bertsekas")
1230e0ed867bSAlp Dener . -as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
1231e0ed867bSAlp Dener . -as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
1232e0ed867bSAlp Dener . -sval - (developer) Hessian perturbation starting value
1233e0ed867bSAlp Dener . -imin - (developer) minimum initial Hessian perturbation
1234e0ed867bSAlp Dener . -imax - (developer) maximum initial Hessian perturbation
1235e0ed867bSAlp Dener . -pmin - (developer) minimum Hessian perturbation
1236e0ed867bSAlp Dener . -pmax - (developer) aximum Hessian perturbation
1237e0ed867bSAlp Dener . -pgfac - (developer) Hessian perturbation growth factor
1238e0ed867bSAlp Dener . -psfac - (developer) Hessian perturbation shrink factor
1239e0ed867bSAlp Dener . -imfac - (developer) initial merit factor for Hessian perturbation
1240e0ed867bSAlp Dener . -pmgfac - (developer) merit growth factor for Hessian perturbation
1241e0ed867bSAlp Dener . -pmsfac - (developer) merit shrink factor for Hessian perturbation
1242e0ed867bSAlp Dener . -eta1 - (developer) threshold for rejecting step (-update_type reduction)
1243e0ed867bSAlp Dener . -eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
1244e0ed867bSAlp Dener . -eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
1245e0ed867bSAlp Dener . -eta4 - (developer) threshold for accepting good step (-update_type reduction)
1246e0ed867bSAlp Dener . -alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
1247e0ed867bSAlp Dener . -alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
1248e0ed867bSAlp Dener . -alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
1249e0ed867bSAlp Dener . -alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
1250e0ed867bSAlp Dener . -alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
1251e0ed867bSAlp Dener . -epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
1252e0ed867bSAlp Dener . -mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
1253e0ed867bSAlp Dener . -mu2 - (developer) threshold for accepting good step (-update_type interpolation)
1254e0ed867bSAlp Dener . -gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
1255e0ed867bSAlp Dener . -gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
1256e0ed867bSAlp Dener . -gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
1257e0ed867bSAlp Dener . -gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
1258e0ed867bSAlp Dener . -theta - (developer) trust region interpolation factor (-update_type interpolation)
1259e0ed867bSAlp Dener . -nu1 - (developer) threshold for small line-search step length (-update_type step)
1260e0ed867bSAlp Dener . -nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
1261e0ed867bSAlp Dener . -nu3 - (developer) threshold for large line-search step length (-update_type step)
1262e0ed867bSAlp Dener . -nu4 - (developer) threshold for very large line-search step length (-update_type step)
1263e0ed867bSAlp Dener . -omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
1264e0ed867bSAlp Dener . -omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
1265e0ed867bSAlp Dener . -omega3 - (developer) radius factor for decent line-search step length (-update_type step)
1266e0ed867bSAlp Dener . -omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
1267e0ed867bSAlp Dener . -omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
1268e0ed867bSAlp Dener . -mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
1269e0ed867bSAlp Dener . -mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
1270e0ed867bSAlp Dener . -gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
1271e0ed867bSAlp Dener . -gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
1272e0ed867bSAlp Dener . -gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
1273e0ed867bSAlp Dener . -gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
1274e0ed867bSAlp Dener - -theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1275eb910715SAlp Dener 
1276eb910715SAlp Dener   Level: beginner
1277eb910715SAlp Dener M*/
1278eb910715SAlp Dener 
1279eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1280eb910715SAlp Dener {
1281eb910715SAlp Dener   TAO_BNK        *bnk;
1282eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1283eb910715SAlp Dener   PetscErrorCode ierr;
1284b9ac7092SAlp Dener   PC             pc;
1285eb910715SAlp Dener 
1286eb910715SAlp Dener   PetscFunctionBegin;
1287eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1288eb910715SAlp Dener 
1289eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1290eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1291eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1292eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1293eb910715SAlp Dener 
1294eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1295eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1296eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1297eb910715SAlp Dener 
1298eb910715SAlp Dener   tao->data = (void*)bnk;
1299eb910715SAlp Dener 
130066ed3702SAlp Dener   /*  Hessian shifting parameters */
1301e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1302e0ed867bSAlp Dener   bnk->computestep = TaoBNKComputeStep;
1303e0ed867bSAlp Dener 
1304eb910715SAlp Dener   bnk->sval   = 0.0;
1305eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1306eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1307eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1308eb910715SAlp Dener 
1309eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1310eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1311eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1312eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1313eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1314eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1315eb910715SAlp Dener 
1316eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1317eb910715SAlp Dener   bnk->nu1 = 0.25;
1318eb910715SAlp Dener   bnk->nu2 = 0.50;
1319eb910715SAlp Dener   bnk->nu3 = 1.00;
1320eb910715SAlp Dener   bnk->nu4 = 1.25;
1321eb910715SAlp Dener 
1322eb910715SAlp Dener   bnk->omega1 = 0.25;
1323eb910715SAlp Dener   bnk->omega2 = 0.50;
1324eb910715SAlp Dener   bnk->omega3 = 1.00;
1325eb910715SAlp Dener   bnk->omega4 = 2.00;
1326eb910715SAlp Dener   bnk->omega5 = 4.00;
1327eb910715SAlp Dener 
1328eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1329eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1330eb910715SAlp Dener   bnk->eta2 = 0.25;
1331eb910715SAlp Dener   bnk->eta3 = 0.50;
1332eb910715SAlp Dener   bnk->eta4 = 0.90;
1333eb910715SAlp Dener 
1334eb910715SAlp Dener   bnk->alpha1 = 0.25;
1335eb910715SAlp Dener   bnk->alpha2 = 0.50;
1336eb910715SAlp Dener   bnk->alpha3 = 1.00;
1337eb910715SAlp Dener   bnk->alpha4 = 2.00;
1338eb910715SAlp Dener   bnk->alpha5 = 4.00;
1339eb910715SAlp Dener 
1340eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1341eb910715SAlp Dener   bnk->mu1 = 0.10;
1342eb910715SAlp Dener   bnk->mu2 = 0.50;
1343eb910715SAlp Dener 
1344eb910715SAlp Dener   bnk->gamma1 = 0.25;
1345eb910715SAlp Dener   bnk->gamma2 = 0.50;
1346eb910715SAlp Dener   bnk->gamma3 = 2.00;
1347eb910715SAlp Dener   bnk->gamma4 = 4.00;
1348eb910715SAlp Dener 
1349eb910715SAlp Dener   bnk->theta = 0.05;
1350eb910715SAlp Dener 
1351eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1352eb910715SAlp Dener   bnk->mu1_i = 0.35;
1353eb910715SAlp Dener   bnk->mu2_i = 0.50;
1354eb910715SAlp Dener 
1355eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1356eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1357eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1358eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1359eb910715SAlp Dener 
1360eb910715SAlp Dener   bnk->theta_i = 0.25;
1361eb910715SAlp Dener 
1362eb910715SAlp Dener   /*  Remaining parameters */
1363c0f10754SAlp Dener   bnk->max_cg_its = 0;
1364eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1365eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1366770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
13670a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
13680a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
136962675beeSAlp Dener   bnk->dmin = 1.0e-6;
137062675beeSAlp Dener   bnk->dmax = 1.0e6;
1371eb910715SAlp Dener 
1372b9ac7092SAlp Dener   bnk->M               = 0;
1373b9ac7092SAlp Dener   bnk->bfgs_pre        = 0;
1374eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
13757b1c7716SAlp Dener   bnk->update_type     = BNK_UPDATE_REDUCTION;
13762f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1377eb910715SAlp Dener 
1378e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1379e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1380e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1381e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1382e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1383e031d6f5SAlp Dener 
1384c0f10754SAlp Dener   /* Create the line search */
1385eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1386eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1387e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1388eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1389eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1390eb910715SAlp Dener 
1391eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1392eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1393eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1394e0ed867bSAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,"tao_bnk_");CHKERRQ(ierr);
1395eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1396b9ac7092SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);
1397b9ac7092SAlp Dener   ierr = PCSetType(pc, PCLMVM);CHKERRQ(ierr);
1398eb910715SAlp Dener   PetscFunctionReturn(0);
1399eb910715SAlp Dener }
1400