xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener #include <petscksp.h>
4eb910715SAlp Dener 
570a3f44bSAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
670a3f44bSAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
770a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
870a3f44bSAlp Dener 
9e031d6f5SAlp Dener /*------------------------------------------------------------*/
10e031d6f5SAlp Dener 
11df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
12df278d8fSAlp Dener 
13c0f10754SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
14eb910715SAlp Dener {
15eb910715SAlp Dener   PetscErrorCode    ierr;
16eb910715SAlp Dener   TAO_BNK           *bnk = (TAO_BNK *)tao->data;
17eb910715SAlp Dener   PC                pc;
1889da521bSAlp Dener   PetscReal         f_min, ftrial, prered, actred, kappa, sigma, resnorm;
19eb910715SAlp Dener   PetscReal         tau, tau_1, tau_2, tau_max, tau_min, max_radius;
200ad3a497SAlp Dener   PetscBool         is_bfgs, is_jacobi, is_symmetric, sym_set;
21c4b75bccSAlp Dener   PetscInt          n, N, nDiff;
22eb910715SAlp Dener   PetscInt          i_max = 5;
23eb910715SAlp Dener   PetscInt          j_max = 1;
24eb910715SAlp Dener   PetscInt          i, j;
252e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
26eb910715SAlp Dener 
27eb910715SAlp Dener   PetscFunctionBegin;
2828017e9fSAlp Dener   /* Project the current point onto the feasible set */
2928017e9fSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
30e031d6f5SAlp Dener   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
31b9ac7092SAlp Dener   if (tao->bounded) {
3228017e9fSAlp Dener     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
33cd929ea3SAlp Dener   }
3428017e9fSAlp Dener 
3528017e9fSAlp Dener   /* Project the initial point onto the feasible region */
363b063059SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
3728017e9fSAlp Dener 
3828017e9fSAlp Dener   /* Check convergence criteria */
3928017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
4061be54a6SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
4161be54a6SAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
4261be54a6SAlp Dener   ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
43f5766c09SAlp Dener   ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
4428017e9fSAlp Dener 
45c0f10754SAlp Dener   /* Test the initial point for convergence */
4689da521bSAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
4789da521bSAlp Dener   ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
48691b26d3SBarry Smith   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
49e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
50e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
51c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
52c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
53c0f10754SAlp Dener 
54e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
55eb910715SAlp Dener   bnk->ksp_atol = 0;
56eb910715SAlp Dener   bnk->ksp_rtol = 0;
57eb910715SAlp Dener   bnk->ksp_dtol = 0;
58eb910715SAlp Dener   bnk->ksp_ctol = 0;
59eb910715SAlp Dener   bnk->ksp_negc = 0;
60eb910715SAlp Dener   bnk->ksp_iter = 0;
61eb910715SAlp Dener   bnk->ksp_othr = 0;
62eb910715SAlp Dener 
63e031d6f5SAlp Dener   /* Reset accepted step type counters */
64e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
65e031d6f5SAlp Dener   bnk->newt = 0;
66e031d6f5SAlp Dener   bnk->bfgs = 0;
67e031d6f5SAlp Dener   bnk->sgrad = 0;
68e031d6f5SAlp Dener   bnk->grad = 0;
69e031d6f5SAlp Dener 
70fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
71fed79b8eSAlp Dener   bnk->pert = bnk->sval;
72fed79b8eSAlp Dener 
73937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
74937a31a1SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
75937a31a1SAlp Dener 
76e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
77b9ac7092SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
78b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr);
79b9ac7092SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr);
80b9ac7092SAlp Dener   if (is_bfgs) {
81b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
82b9ac7092SAlp Dener     ierr = PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M);CHKERRQ(ierr);
83eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr);
84eb910715SAlp Dener     ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr);
85b9ac7092SAlp Dener     ierr = MatSetSizes(bnk->M, n, n, N, N);CHKERRQ(ierr);
86cd929ea3SAlp Dener     ierr = MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
870ad3a497SAlp Dener     ierr = MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric);CHKERRQ(ierr);
880ad3a497SAlp Dener     if (!sym_set || !is_symmetric) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
89b9ac7092SAlp Dener   } else if (is_jacobi) {
90b9ac7092SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
91e031d6f5SAlp Dener   }
92e031d6f5SAlp Dener 
93e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
9462675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
9562675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
96eb910715SAlp Dener 
97eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
98eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
99c0f10754SAlp Dener   *needH = PETSC_TRUE;
1002e6e4ca1SStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)tao->ksp,"KSPCGSetRadius_C",&kspTR);CHKERRQ(ierr);
1012e6e4ca1SStefano Zampini   if (kspTR) {
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 */
1189fa2b5dcSStefano Zampini           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 {
124c5e9d94cSAlp Dener             ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
125c5e9d94cSAlp Dener             bnk->H_inactive = tao->hessian;
12628017e9fSAlp Dener           }
127c0f10754SAlp Dener           *needH = PETSC_FALSE;
128eb910715SAlp Dener         }
129eb910715SAlp Dener 
130eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
13162602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
13262602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
13362602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
1343b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
13589da521bSAlp Dener           /* Compute the step we actually accepted */
136eb910715SAlp Dener           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
13762602cfbSAlp Dener           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
13862602cfbSAlp Dener           /* Compute the objective at the trial */
13962602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
140691b26d3SBarry Smith           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
14162602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
142eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
143eb910715SAlp Dener             tau = bnk->gamma1_i;
144eb910715SAlp Dener           } else {
1450a4511e9SAlp Dener             if (ftrial < f_min) {
1460a4511e9SAlp Dener               f_min = ftrial;
147eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
148eb910715SAlp Dener             }
14908752603SAlp Dener 
150770b7498SAlp Dener             /* Compute the predicted and actual reduction */
15189da521bSAlp Dener             if (bnk->active_idx) {
15208752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
15308752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1542ab2a32cSAlp Dener             } else {
15508752603SAlp Dener               bnk->X_inactive = bnk->W;
15608752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1572ab2a32cSAlp Dener             }
15808752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
15908752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
16089da521bSAlp Dener             if (bnk->active_idx) {
16108752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
16208752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1632ab2a32cSAlp Dener             }
164eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
165eb910715SAlp Dener             actred = bnk->f - ftrial;
1663105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
167eb910715SAlp Dener               kappa = 1.0;
1683105154fSTodd Munson             } else {
169eb910715SAlp Dener               kappa = actred / prered;
170eb910715SAlp Dener             }
171eb910715SAlp Dener 
172eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
173eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
174eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
175eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
176eb910715SAlp Dener 
17718cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
178eb910715SAlp Dener               /*  Great agreement */
179eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
180eb910715SAlp Dener 
181eb910715SAlp Dener               if (tau_max < 1.0) {
182eb910715SAlp Dener                 tau = bnk->gamma3_i;
1833105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
184eb910715SAlp Dener                 tau = bnk->gamma4_i;
1853105154fSTodd Munson               } else {
186eb910715SAlp Dener                 tau = tau_max;
187eb910715SAlp Dener               }
18818cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
189eb910715SAlp Dener               /*  Good agreement */
190eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
191eb910715SAlp Dener 
192eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
193eb910715SAlp Dener                 tau = bnk->gamma2_i;
194eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
195eb910715SAlp Dener                 tau = bnk->gamma3_i;
196eb910715SAlp Dener               } else {
197eb910715SAlp Dener                 tau = tau_max;
198eb910715SAlp Dener               }
1998f8a4e06SAlp Dener             } else {
200eb910715SAlp Dener               /*  Not good agreement */
201eb910715SAlp Dener               if (tau_min > 1.0) {
202eb910715SAlp Dener                 tau = bnk->gamma2_i;
203eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
204eb910715SAlp Dener                 tau = bnk->gamma1_i;
205eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
206eb910715SAlp Dener                 tau = bnk->gamma1_i;
2073105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
208eb910715SAlp Dener                 tau = tau_1;
2093105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
210eb910715SAlp Dener                 tau = tau_2;
211eb910715SAlp Dener               } else {
212eb910715SAlp Dener                 tau = tau_max;
213eb910715SAlp Dener               }
214eb910715SAlp Dener             }
215eb910715SAlp Dener           }
216eb910715SAlp Dener           tao->trust = tau * tao->trust;
217eb910715SAlp Dener         }
218eb910715SAlp Dener 
2190a4511e9SAlp Dener         if (f_min < bnk->f) {
220937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2210a4511e9SAlp Dener           bnk->f = f_min;
222937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
223eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
2243b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
225937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
226937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
22709164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
22861be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
22961be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
23061be54a6SAlp Dener           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
231937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
232f5766c09SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
233c0f10754SAlp Dener           *needH = PETSC_TRUE;
234937a31a1SAlp Dener           /* Test the new step for convergence */
23589da521bSAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
23689da521bSAlp Dener           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
237691b26d3SBarry Smith           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
238e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
239e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
240eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
241eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
242937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
243414d97d3SAlp Dener           ierr = TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
244eb910715SAlp Dener         }
245eb910715SAlp Dener       }
246eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
247e031d6f5SAlp Dener 
248e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
249e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
250e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
251eb910715SAlp Dener       break;
252eb910715SAlp Dener 
253eb910715SAlp Dener     default:
254eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
255eb910715SAlp Dener       tao->trust = 0.0;
256eb910715SAlp Dener       break;
257eb910715SAlp Dener     }
258eb910715SAlp Dener   }
259eb910715SAlp Dener   PetscFunctionReturn(0);
260eb910715SAlp Dener }
261eb910715SAlp Dener 
262df278d8fSAlp Dener /*------------------------------------------------------------*/
263df278d8fSAlp Dener 
264e0ed867bSAlp Dener /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */
26562675beeSAlp Dener 
26662675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
26762675beeSAlp Dener {
26862675beeSAlp Dener   PetscErrorCode ierr;
26962675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
27062675beeSAlp Dener 
27162675beeSAlp Dener   PetscFunctionBegin;
27262675beeSAlp Dener   /* Compute the Hessian */
27362675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
27462675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
275b9ac7092SAlp Dener   if (bnk->M) {
27662675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
27762675beeSAlp Dener   }
278e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
279f5766c09SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
280e0ed867bSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
281f5766c09SAlp Dener   if (bnk->active_idx) {
282e0ed867bSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
283e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
284f5766c09SAlp Dener       ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
285e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
286e0ed867bSAlp Dener     } else {
287e0ed867bSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
288e0ed867bSAlp Dener     }
289e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
290e0ed867bSAlp Dener       ierr = PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx);CHKERRQ(ierr);
291e0ed867bSAlp Dener     }
292e0ed867bSAlp Dener   } else {
293c5e9d94cSAlp Dener     ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
294c5e9d94cSAlp Dener     bnk->H_inactive = tao->hessian;
295e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
296f5766c09SAlp Dener       ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
297e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
298e0ed867bSAlp Dener     } else {
299c5e9d94cSAlp Dener       ierr = PetscObjectReference((PetscObject)tao->hessian_pre);
300c5e9d94cSAlp Dener       bnk->Hpre_inactive = tao->hessian_pre;
301e0ed867bSAlp Dener     }
302e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
303e0ed867bSAlp Dener       ierr = PCLMVMClearIS(bnk->bfgs_pre);CHKERRQ(ierr);
304e0ed867bSAlp Dener     }
305e0ed867bSAlp Dener   }
30662675beeSAlp Dener   PetscFunctionReturn(0);
30762675beeSAlp Dener }
30862675beeSAlp Dener 
30962675beeSAlp Dener /*------------------------------------------------------------*/
31062675beeSAlp Dener 
3112f75a4aaSAlp Dener /* Routine for estimating the active set */
3122f75a4aaSAlp Dener 
31308752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3142f75a4aaSAlp Dener {
3152f75a4aaSAlp Dener   PetscErrorCode ierr;
3162f75a4aaSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
317c4b75bccSAlp Dener   PetscBool      hessComputed, diagExists;
3182f75a4aaSAlp Dener 
3192f75a4aaSAlp Dener   PetscFunctionBegin;
32008752603SAlp Dener   switch (asType) {
3212f75a4aaSAlp Dener   case BNK_AS_NONE:
3222f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3232f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3242f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3252f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3262f75a4aaSAlp Dener     break;
3272f75a4aaSAlp Dener 
3282f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3292f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
330b9ac7092SAlp Dener     if (bnk->M) {
3312f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3329515a401SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3332f75a4aaSAlp Dener     } else {
334fc5ca067SStefano Zampini       hessComputed = diagExists = PETSC_FALSE;
335f5766c09SAlp Dener       if (tao->hessian) {
33661be54a6SAlp Dener         ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
337f5766c09SAlp Dener       }
338fc5ca067SStefano Zampini       if (hessComputed) {
339fc5ca067SStefano Zampini         ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
340fc5ca067SStefano Zampini       }
341fc5ca067SStefano Zampini       if (diagExists) {
3429b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3432f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3449b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3459b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
346c3b366b1Sprj-         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);
3472f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
34861be54a6SAlp Dener       } else {
349c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
35061be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
35161be54a6SAlp Dener       }
3522f75a4aaSAlp Dener     }
3532f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
35489da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
35589da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
356c4b75bccSAlp Dener     break;
3572f75a4aaSAlp Dener 
3582f75a4aaSAlp Dener   default:
3592f75a4aaSAlp Dener     break;
3602f75a4aaSAlp Dener   }
3612f75a4aaSAlp Dener   PetscFunctionReturn(0);
3622f75a4aaSAlp Dener }
3632f75a4aaSAlp Dener 
3642f75a4aaSAlp Dener /*------------------------------------------------------------*/
3652f75a4aaSAlp Dener 
3662f75a4aaSAlp Dener /* Routine for bounding the step direction */
3672f75a4aaSAlp Dener 
368a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
3692f75a4aaSAlp Dener {
3702f75a4aaSAlp Dener   PetscErrorCode               ierr;
3712f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3722f75a4aaSAlp Dener 
3732f75a4aaSAlp Dener   PetscFunctionBegin;
374a1318120SAlp Dener   switch (asType) {
3752f75a4aaSAlp Dener   case BNK_AS_NONE:
376c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
3772f75a4aaSAlp Dener     break;
3782f75a4aaSAlp Dener 
3792f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
380c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
3812f75a4aaSAlp Dener     break;
3822f75a4aaSAlp Dener 
3832f75a4aaSAlp Dener   default:
3842f75a4aaSAlp Dener     break;
3852f75a4aaSAlp Dener   }
3862f75a4aaSAlp Dener   PetscFunctionReturn(0);
3872f75a4aaSAlp Dener }
3882f75a4aaSAlp Dener 
389e031d6f5SAlp Dener /*------------------------------------------------------------*/
390e031d6f5SAlp Dener 
391e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
392e031d6f5SAlp Dener    accelerate Newton convergence.
393e031d6f5SAlp Dener 
394e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
395e031d6f5SAlp Dener    for more gradient evaluations.
396e031d6f5SAlp Dener */
397e031d6f5SAlp Dener 
398c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
399c0f10754SAlp Dener {
400c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
401c0f10754SAlp Dener   PetscErrorCode               ierr;
402c0f10754SAlp Dener 
403c0f10754SAlp Dener   PetscFunctionBegin;
404c0f10754SAlp Dener   *terminate = PETSC_FALSE;
405c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
406c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
407c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
408c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
409c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
410c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
411c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
412c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
413c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
414c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
415e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
416c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
417c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
418c0f10754SAlp Dener     if (bnk->bncg->reason == TAO_CONVERGED_GATOL || bnk->bncg->reason == TAO_CONVERGED_GRTOL || bnk->bncg->reason == TAO_CONVERGED_GTTOL || bnk->bncg->reason == TAO_CONVERGED_MINF) {
419c0f10754SAlp Dener       *terminate = PETSC_TRUE;
42061be54a6SAlp Dener     } else {
42133c78596SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
422c0f10754SAlp Dener     }
423c0f10754SAlp Dener   }
424c0f10754SAlp Dener   PetscFunctionReturn(0);
425c0f10754SAlp Dener }
426c0f10754SAlp Dener 
4272f75a4aaSAlp Dener /*------------------------------------------------------------*/
4282f75a4aaSAlp Dener 
429c0f10754SAlp Dener /* Routine for computing the Newton step. */
430df278d8fSAlp Dener 
4316b591159SAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
432eb910715SAlp Dener {
433eb910715SAlp Dener   PetscErrorCode    ierr;
434eb910715SAlp Dener   TAO_BNK           *bnk = (TAO_BNK *)tao->data;
435eb910715SAlp Dener   PetscInt          bfgsUpdates = 0;
436eb910715SAlp Dener   PetscInt          kspits;
437bddd1ffdSAlp Dener   PetscBool         is_lmvm;
4382e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
439eb910715SAlp Dener 
440eb910715SAlp Dener   PetscFunctionBegin;
44189da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
44289da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
44389da521bSAlp Dener   if (!bnk->inactive_idx) {
44489da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
445a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
44689da521bSAlp Dener     PetscFunctionReturn(0);
44789da521bSAlp Dener   }
44889da521bSAlp Dener 
44962675beeSAlp Dener   /* Shift the reduced Hessian matrix */
450e831869dSStefano Zampini   if (shift && bnk->pert > 0) {
451f7bf01afSAlp Dener     ierr = PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm);CHKERRQ(ierr);
452f7bf01afSAlp Dener     if (is_lmvm) {
453f7bf01afSAlp Dener       ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
454f7bf01afSAlp Dener     } else {
45562675beeSAlp Dener       ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
45662675beeSAlp Dener       if (bnk->H_inactive != bnk->Hpre_inactive) {
45762675beeSAlp Dener         ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
45862675beeSAlp Dener       }
45962675beeSAlp Dener     }
460f7bf01afSAlp Dener   }
46162675beeSAlp Dener 
462eb910715SAlp Dener   /* Solve the Newton system of equations */
463937a31a1SAlp Dener   tao->ksp_its = 0;
4642f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
4655e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
46694d5fdc2STristan Konolige   ierr = KSPResetFromOptions(tao->ksp);CHKERRQ(ierr);
46709164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
4685e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
46989da521bSAlp Dener   if (bnk->active_idx) {
4705e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
4715e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
4725e9b73cbSAlp Dener   } else {
4735e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4745e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
47528017e9fSAlp Dener   }
4762e6e4ca1SStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)tao->ksp,"KSPCGSetRadius_C",&kspTR);CHKERRQ(ierr);
4772e6e4ca1SStefano Zampini   if (kspTR) {
478fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
4795e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
480eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
481eb910715SAlp Dener     tao->ksp_its += kspits;
482eb910715SAlp Dener     tao->ksp_tot_its += kspits;
483080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
484eb910715SAlp Dener 
485eb910715SAlp Dener     if (0.0 == tao->trust) {
486eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
487080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
488080d2917SAlp Dener         tao->trust = bnk->dnorm;
489eb910715SAlp Dener 
490eb910715SAlp Dener         /* Modify the radius if it is too large or small */
491eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
492eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
493eb910715SAlp Dener       } else {
494eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
495eb910715SAlp Dener            the trust-region subproblem to get a direction */
496eb910715SAlp Dener         tao->trust = tao->trust0;
497eb910715SAlp Dener 
498eb910715SAlp Dener         /* Modify the radius if it is too large or small */
499eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
500eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
501eb910715SAlp Dener 
502fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5035e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
504eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
505eb910715SAlp Dener         tao->ksp_its += kspits;
506eb910715SAlp Dener         tao->ksp_tot_its += kspits;
507080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
508eb910715SAlp Dener 
509691b26d3SBarry Smith         if (bnk->dnorm == 0.0) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero");
510eb910715SAlp Dener       }
511eb910715SAlp Dener     }
512eb910715SAlp Dener   } else {
5135e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
514eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
515eb910715SAlp Dener     tao->ksp_its += kspits;
516eb910715SAlp Dener     tao->ksp_tot_its += kspits;
517eb910715SAlp Dener   }
5185e9b73cbSAlp Dener   /* Restore sub vectors back */
51989da521bSAlp Dener   if (bnk->active_idx) {
5205e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5215e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5225e9b73cbSAlp Dener   }
523770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
524fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
525a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
526770b7498SAlp Dener 
527770b7498SAlp Dener   /* Record convergence reasons */
528e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
529e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
530770b7498SAlp Dener     ++bnk->ksp_atol;
531e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
532770b7498SAlp Dener     ++bnk->ksp_rtol;
533e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
534770b7498SAlp Dener     ++bnk->ksp_ctol;
535e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
536770b7498SAlp Dener     ++bnk->ksp_negc;
537e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
538770b7498SAlp Dener     ++bnk->ksp_dtol;
539e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
540770b7498SAlp Dener     ++bnk->ksp_iter;
541770b7498SAlp Dener   } else {
542770b7498SAlp Dener     ++bnk->ksp_othr;
543770b7498SAlp Dener   }
544fed79b8eSAlp Dener 
545fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
546b9ac7092SAlp Dener   if (bnk->M) {
547cd929ea3SAlp Dener     ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
548b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
549fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
550cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
55109164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
552eb910715SAlp Dener     }
553fed79b8eSAlp Dener   }
5546b591159SAlp Dener   *step_type = BNK_NEWTON;
555e465cd6fSAlp Dener   PetscFunctionReturn(0);
556e465cd6fSAlp Dener }
557eb910715SAlp Dener 
55862675beeSAlp Dener /*------------------------------------------------------------*/
55962675beeSAlp Dener 
5605e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5615e9b73cbSAlp Dener 
5625e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
5635e9b73cbSAlp Dener {
5645e9b73cbSAlp Dener   PetscErrorCode ierr;
5655e9b73cbSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
5665e9b73cbSAlp Dener 
5675e9b73cbSAlp Dener   PetscFunctionBegin;
5685e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
56989da521bSAlp Dener   if (bnk->active_idx) {
5705e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5715e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5725e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5735e9b73cbSAlp Dener   } else {
5745e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
5755e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5765e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
5775e9b73cbSAlp Dener   }
5785e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5795e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
5805e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
58133c78596SAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);CHKERRQ(ierr);
5825e9b73cbSAlp Dener   /* Restore the sub vectors */
58389da521bSAlp Dener   if (bnk->active_idx) {
5845e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5855e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
5865e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5875e9b73cbSAlp Dener   }
5885e9b73cbSAlp Dener   PetscFunctionReturn(0);
5895e9b73cbSAlp Dener }
5905e9b73cbSAlp Dener 
5915e9b73cbSAlp Dener /*------------------------------------------------------------*/
5925e9b73cbSAlp Dener 
59362675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
59462675beeSAlp Dener 
59562675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
59662675beeSAlp Dener    in the event that the Newton step fails the test.
59762675beeSAlp Dener */
59862675beeSAlp Dener 
599e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
600e465cd6fSAlp Dener {
601e465cd6fSAlp Dener   PetscErrorCode ierr;
602e465cd6fSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
603b2d8c577SAlp Dener   PetscReal      gdx, e_min;
604e465cd6fSAlp Dener   PetscInt       bfgsUpdates;
605e465cd6fSAlp Dener 
606e465cd6fSAlp Dener   PetscFunctionBegin;
6076b591159SAlp Dener   switch (*stepType) {
6086b591159SAlp Dener   case BNK_NEWTON:
609080d2917SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
610eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
611eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
612eb910715SAlp Dener         Update the perturbation for next time */
613eb910715SAlp Dener       if (bnk->pert <= 0.0) {
6142e6e4ca1SStefano Zampini         PetscBool is_gltr;
6152e6e4ca1SStefano Zampini 
616eb910715SAlp Dener         /* Initialize the perturbation */
617eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
6182e6e4ca1SStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr);CHKERRQ(ierr);
6192e6e4ca1SStefano Zampini         if (is_gltr) {
62005de396fSBarry Smith           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
621eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
622eb910715SAlp Dener         }
623eb910715SAlp Dener       } else {
624eb910715SAlp Dener         /* Increase the perturbation */
625eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
626eb910715SAlp Dener       }
627eb910715SAlp Dener 
6280ad3a497SAlp Dener       if (!bnk->M) {
629eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
630eb910715SAlp Dener           Must use gradient direction in this case */
631080d2917SAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
632eb910715SAlp Dener         *stepType = BNK_GRADIENT;
633eb910715SAlp Dener       } else {
634eb910715SAlp Dener         /* Attempt to use the BFGS direction */
6359515a401SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
636eb910715SAlp Dener 
6378d5ead36SAlp Dener         /* Check for success (descent direction)
6388d5ead36SAlp Dener           NOTE: Negative gdx here means not a descent direction because
6398d5ead36SAlp Dener           the fall-back step is missing a negative sign. */
640080d2917SAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
6413105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
642eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
643eb910715SAlp Dener             We can assert bfgsUpdates > 1 in this case because
644eb910715SAlp Dener             the first solve produces the scaled gradient direction,
645eb910715SAlp Dener             which is guaranteed to be descent */
646eb910715SAlp Dener 
647eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
648cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
64909164190SAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
6509515a401SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
651eb910715SAlp Dener 
652eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
653eb910715SAlp Dener         } else {
654cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
655eb910715SAlp Dener           if (1 == bfgsUpdates) {
656eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
657eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
658eb910715SAlp Dener           } else {
659eb910715SAlp Dener             *stepType = BNK_BFGS;
660eb910715SAlp Dener           }
661eb910715SAlp Dener         }
662eb910715SAlp Dener       }
6638d5ead36SAlp Dener       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6648d5ead36SAlp Dener       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
665a1318120SAlp Dener       ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
666eb910715SAlp Dener     } else {
667eb910715SAlp Dener       /* Computed Newton step is descent */
668eb910715SAlp Dener       switch (ksp_reason) {
669eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
670eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
671eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
672eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
673eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
674eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
675eb910715SAlp Dener         if (bnk->pert <= 0.0) {
6762e6e4ca1SStefano Zampini           PetscBool is_gltr;
6772e6e4ca1SStefano Zampini 
678eb910715SAlp Dener           /* Initialize the perturbation */
679eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
6802e6e4ca1SStefano Zampini           ierr = PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr);CHKERRQ(ierr);
6812e6e4ca1SStefano Zampini           if (is_gltr) {
68205de396fSBarry Smith             ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
683eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
684eb910715SAlp Dener           }
685eb910715SAlp Dener         } else {
686eb910715SAlp Dener           /* Increase the perturbation */
687eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
688eb910715SAlp Dener         }
689eb910715SAlp Dener         break;
690eb910715SAlp Dener 
691eb910715SAlp Dener       default:
692eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
693eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
694eb910715SAlp Dener         if (bnk->pert < bnk->pmin) {
695eb910715SAlp Dener           bnk->pert = 0.0;
696eb910715SAlp Dener         }
697eb910715SAlp Dener         break;
698eb910715SAlp Dener       }
699fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
700eb910715SAlp Dener     }
7016b591159SAlp Dener     break;
7026b591159SAlp Dener 
7036b591159SAlp Dener   case BNK_BFGS:
7046b591159SAlp Dener     /* Check for success (descent direction) */
7056b591159SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
7066b591159SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
7076b591159SAlp Dener       /* Step is not descent or solve was not successful
7086b591159SAlp Dener          Use steepest descent direction (scaled) */
7096b591159SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
7106b591159SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
7119515a401SAlp Dener       ierr = MatSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
7126b591159SAlp Dener       ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
7136b591159SAlp Dener       ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
7146b591159SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
7156b591159SAlp Dener     } else {
7166b591159SAlp Dener       *stepType = BNK_BFGS;
7176b591159SAlp Dener     }
7186b591159SAlp Dener     break;
7196b591159SAlp Dener 
7206b591159SAlp Dener   case BNK_SCALED_GRADIENT:
7216b591159SAlp Dener     break;
7226b591159SAlp Dener 
7236b591159SAlp Dener   default:
7246b591159SAlp Dener     break;
7256b591159SAlp Dener   }
7266b591159SAlp Dener 
727eb910715SAlp Dener   PetscFunctionReturn(0);
728eb910715SAlp Dener }
729eb910715SAlp Dener 
730df278d8fSAlp Dener /*------------------------------------------------------------*/
731df278d8fSAlp Dener 
732df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
733df278d8fSAlp Dener 
734df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
735df278d8fSAlp Dener   Newton step does not produce a valid step length.
736df278d8fSAlp Dener */
737df278d8fSAlp Dener 
738937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
739c14b763aSAlp Dener {
740c14b763aSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
741c14b763aSAlp Dener   PetscErrorCode               ierr;
742c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
743b2d8c577SAlp Dener   PetscReal                    e_min, gdx;
744c14b763aSAlp Dener   PetscInt                     bfgsUpdates;
745c14b763aSAlp Dener 
746c14b763aSAlp Dener   PetscFunctionBegin;
747c14b763aSAlp Dener   /* Perform the linesearch */
748c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
749c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
750c14b763aSAlp Dener 
751b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
752c14b763aSAlp Dener     /* Linesearch failed, revert solution */
753c14b763aSAlp Dener     bnk->f = bnk->fold;
754c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
755c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
756c14b763aSAlp Dener 
757937a31a1SAlp Dener     switch(*stepType) {
758c14b763aSAlp Dener     case BNK_NEWTON:
7598d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
760c14b763aSAlp Dener          Update the perturbation for next time */
761c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
7622e6e4ca1SStefano Zampini         PetscBool is_gltr;
7632e6e4ca1SStefano Zampini 
764c14b763aSAlp Dener         /* Initialize the perturbation */
765c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
7662e6e4ca1SStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr);CHKERRQ(ierr);
7672e6e4ca1SStefano Zampini         if (is_gltr) {
76805de396fSBarry Smith           ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
769c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
770c14b763aSAlp Dener         }
771c14b763aSAlp Dener       } else {
772c14b763aSAlp Dener         /* Increase the perturbation */
773c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
774c14b763aSAlp Dener       }
775c14b763aSAlp Dener 
7760ad3a497SAlp Dener       if (!bnk->M) {
777c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
778c14b763aSAlp Dener            Must use gradient direction in this case */
779937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
780937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
781c14b763aSAlp Dener       } else {
782c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
7839515a401SAlp Dener         ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
7848d5ead36SAlp Dener         /* Check for success (descent direction)
7858d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
786c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7873105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
788c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
789c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
790c14b763aSAlp Dener              Use steepest descent direction (scaled) */
791cd929ea3SAlp Dener           ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
792c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
7939515a401SAlp Dener           ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
794c14b763aSAlp Dener 
795c14b763aSAlp Dener           bfgsUpdates = 1;
796937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
797c14b763aSAlp Dener         } else {
798cd929ea3SAlp Dener           ierr = MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
799c14b763aSAlp Dener           if (1 == bfgsUpdates) {
800c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
801937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
802c14b763aSAlp Dener           } else {
803937a31a1SAlp Dener             *stepType = BNK_BFGS;
804c14b763aSAlp Dener           }
805c14b763aSAlp Dener         }
806c14b763aSAlp Dener       }
807c14b763aSAlp Dener       break;
808c14b763aSAlp Dener 
809c14b763aSAlp Dener     case BNK_BFGS:
810c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
811c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
812c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
813cd929ea3SAlp Dener       ierr = MatLMVMReset(bnk->M, PETSC_FALSE);CHKERRQ(ierr);
814c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
8159515a401SAlp Dener       ierr = MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
816c14b763aSAlp Dener 
817c14b763aSAlp Dener       bfgsUpdates = 1;
818937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
819c14b763aSAlp Dener       break;
820c14b763aSAlp Dener     }
8218d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8228d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
823a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
824c14b763aSAlp Dener 
8258d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
826c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
827c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
828c14b763aSAlp Dener   }
829c14b763aSAlp Dener   *reason = ls_reason;
830c14b763aSAlp Dener   PetscFunctionReturn(0);
831c14b763aSAlp Dener }
832c14b763aSAlp Dener 
833df278d8fSAlp Dener /*------------------------------------------------------------*/
834df278d8fSAlp Dener 
835df278d8fSAlp Dener /* Routine for updating the trust radius.
836df278d8fSAlp Dener 
837df278d8fSAlp Dener   Function features three different update methods:
838df278d8fSAlp Dener   1) Line-search step length based
839df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
840df278d8fSAlp Dener   3) Interpolation
841df278d8fSAlp Dener */
842df278d8fSAlp Dener 
84328017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
844080d2917SAlp Dener {
845080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
846080d2917SAlp Dener   PetscErrorCode ierr;
847080d2917SAlp Dener 
848b1c2d0e3SAlp Dener   PetscReal      step, kappa;
849080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
850080d2917SAlp Dener 
851080d2917SAlp Dener   PetscFunctionBegin;
852080d2917SAlp Dener   /* Update trust region radius */
853080d2917SAlp Dener   *accept = PETSC_FALSE;
85428017e9fSAlp Dener   switch(updateType) {
855080d2917SAlp Dener   case BNK_UPDATE_STEP:
856c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
857080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
858080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
859080d2917SAlp Dener       if (step < bnk->nu1) {
860080d2917SAlp Dener         /* Very bad step taken; reduce radius */
861080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
862080d2917SAlp Dener       } else if (step < bnk->nu2) {
863080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
864080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
865080d2917SAlp Dener       } else if (step < bnk->nu3) {
866080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
867080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
868080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
869080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
870080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
871080d2917SAlp Dener         }
872080d2917SAlp Dener       } else if (step < bnk->nu4) {
873080d2917SAlp Dener         /*  Full step taken; increase the radius */
874080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
875080d2917SAlp Dener       } else {
876080d2917SAlp Dener         /*  More than full step taken; increase the radius */
877080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
878080d2917SAlp Dener       }
879080d2917SAlp Dener     } else {
880080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
881080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
882080d2917SAlp Dener     }
883080d2917SAlp Dener     break;
884080d2917SAlp Dener 
885080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
886080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
887e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
888fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
889fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
890fed79b8eSAlp Dener            be rejected! */
891080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8923105154fSTodd Munson       } else {
893b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
894080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
895080d2917SAlp Dener         } else {
8963105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
897080d2917SAlp Dener             kappa = 1.0;
8983105154fSTodd Munson           } else {
899080d2917SAlp Dener             kappa = actred / prered;
900080d2917SAlp Dener           }
901fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
902080d2917SAlp Dener           if (kappa < bnk->eta1) {
903fed79b8eSAlp Dener             /* Reject the step */
904080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
9053105154fSTodd Munson           } else {
906fed79b8eSAlp Dener             /* Accept the step */
907c133c014SAlp Dener             *accept = PETSC_TRUE;
908c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
9098d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
910080d2917SAlp Dener               if (kappa < bnk->eta2) {
911080d2917SAlp Dener                 /* Marginal bad step */
912c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
9133105154fSTodd Munson               } else if (kappa < bnk->eta3) {
914fed79b8eSAlp Dener                 /* Reasonable step */
915fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
9163105154fSTodd Munson               } else if (kappa < bnk->eta4) {
917080d2917SAlp Dener                 /* Good step */
918c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
9193105154fSTodd Munson               } else {
920080d2917SAlp Dener                 /* Very good step */
921c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
922080d2917SAlp Dener               }
923c133c014SAlp Dener             }
924080d2917SAlp Dener           }
925080d2917SAlp Dener         }
926080d2917SAlp Dener       }
927080d2917SAlp Dener     } else {
928080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
929080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
930080d2917SAlp Dener     }
931080d2917SAlp Dener     break;
932080d2917SAlp Dener 
933080d2917SAlp Dener   default:
934080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
935b1c2d0e3SAlp Dener       if (prered < 0.0) {
936080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
937080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
938080d2917SAlp Dener         /*  be rejected! */
939080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
940080d2917SAlp Dener       } else {
941b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
942080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
943080d2917SAlp Dener         } else {
944080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
945080d2917SAlp Dener             kappa = 1.0;
946080d2917SAlp Dener           } else {
947080d2917SAlp Dener             kappa = actred / prered;
948080d2917SAlp Dener           }
949080d2917SAlp Dener 
950080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
951080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
952080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
953080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
954080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
955080d2917SAlp Dener 
956080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
957080d2917SAlp Dener             /*  Great agreement */
958080d2917SAlp Dener             *accept = PETSC_TRUE;
959080d2917SAlp Dener             if (tau_max < 1.0) {
960080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
961080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
962080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
963080d2917SAlp Dener             } else {
964080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
965080d2917SAlp Dener             }
966080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
967080d2917SAlp Dener             /*  Good agreement */
968080d2917SAlp Dener             *accept = PETSC_TRUE;
969080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
970080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
971080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
972080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
973080d2917SAlp Dener             } else if (tau_max < 1.0) {
974080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
975080d2917SAlp Dener             } else {
976080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
977080d2917SAlp Dener             }
978080d2917SAlp Dener           } else {
979080d2917SAlp Dener             /*  Not good agreement */
980080d2917SAlp Dener             if (tau_min > 1.0) {
981080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
982080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
983080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
984080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
985080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
986080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
987080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
988080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
989080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
990080d2917SAlp Dener             } else {
991080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
992080d2917SAlp Dener             }
993080d2917SAlp Dener           }
994080d2917SAlp Dener         }
995080d2917SAlp Dener       }
996080d2917SAlp Dener     } else {
997080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
998080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
999080d2917SAlp Dener     }
100028017e9fSAlp Dener     break;
1001080d2917SAlp Dener   }
1002c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
1003c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1004fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1005080d2917SAlp Dener   PetscFunctionReturn(0);
1006080d2917SAlp Dener }
1007080d2917SAlp Dener 
1008eb910715SAlp Dener /* ---------------------------------------------------------- */
1009df278d8fSAlp Dener 
101062675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
101162675beeSAlp Dener {
101262675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
101362675beeSAlp Dener 
101462675beeSAlp Dener   PetscFunctionBegin;
101562675beeSAlp Dener   switch (stepType) {
101662675beeSAlp Dener   case BNK_NEWTON:
101762675beeSAlp Dener     ++bnk->newt;
101862675beeSAlp Dener     break;
101962675beeSAlp Dener   case BNK_BFGS:
102062675beeSAlp Dener     ++bnk->bfgs;
102162675beeSAlp Dener     break;
102262675beeSAlp Dener   case BNK_SCALED_GRADIENT:
102362675beeSAlp Dener     ++bnk->sgrad;
102462675beeSAlp Dener     break;
102562675beeSAlp Dener   case BNK_GRADIENT:
102662675beeSAlp Dener     ++bnk->grad;
102762675beeSAlp Dener     break;
102862675beeSAlp Dener   default:
102962675beeSAlp Dener     break;
103062675beeSAlp Dener   }
103162675beeSAlp Dener   PetscFunctionReturn(0);
103262675beeSAlp Dener }
103362675beeSAlp Dener 
103462675beeSAlp Dener /* ---------------------------------------------------------- */
103562675beeSAlp Dener 
10369b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1037eb910715SAlp Dener {
1038eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1039eb910715SAlp Dener   PetscErrorCode ierr;
1040e031d6f5SAlp Dener   PetscInt       i;
1041eb910715SAlp Dener 
1042eb910715SAlp Dener   PetscFunctionBegin;
1043c4b75bccSAlp Dener   if (!tao->gradient) {
1044c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1045c4b75bccSAlp Dener   }
1046c4b75bccSAlp Dener   if (!tao->stepdirection) {
1047c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1048c4b75bccSAlp Dener   }
1049c4b75bccSAlp Dener   if (!bnk->W) {
1050c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1051c4b75bccSAlp Dener   }
1052c4b75bccSAlp Dener   if (!bnk->Xold) {
1053c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1054c4b75bccSAlp Dener   }
1055c4b75bccSAlp Dener   if (!bnk->Gold) {
1056c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1057c4b75bccSAlp Dener   }
1058c4b75bccSAlp Dener   if (!bnk->Xwork) {
1059c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1060c4b75bccSAlp Dener   }
1061c4b75bccSAlp Dener   if (!bnk->Gwork) {
1062c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1063c4b75bccSAlp Dener   }
1064c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1065c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1066c4b75bccSAlp Dener   }
1067c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1068c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1069c4b75bccSAlp Dener   }
1070c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1071c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1072c4b75bccSAlp Dener   }
1073c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1074c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1075c4b75bccSAlp Dener   }
1076e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1077c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1078c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
107989da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
108089da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
108189da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1082c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1083c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1084c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1085c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1086c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1087c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1088c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1089c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1090c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1091c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1092c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1093c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1094c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1095c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1096e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1097e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1098e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1099937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1100e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1101e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1102e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1103e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1104c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1105e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1106e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1107e031d6f5SAlp Dener     }
1108e031d6f5SAlp Dener   }
110983c8fe1dSLisandro Dalcin   bnk->X_inactive = NULL;
111083c8fe1dSLisandro Dalcin   bnk->G_inactive = NULL;
111183c8fe1dSLisandro Dalcin   bnk->inactive_work = NULL;
111283c8fe1dSLisandro Dalcin   bnk->active_work = NULL;
111383c8fe1dSLisandro Dalcin   bnk->inactive_idx = NULL;
111483c8fe1dSLisandro Dalcin   bnk->active_idx = NULL;
111583c8fe1dSLisandro Dalcin   bnk->active_lower = NULL;
111683c8fe1dSLisandro Dalcin   bnk->active_upper = NULL;
111783c8fe1dSLisandro Dalcin   bnk->active_fixed = NULL;
111883c8fe1dSLisandro Dalcin   bnk->M = NULL;
111983c8fe1dSLisandro Dalcin   bnk->H_inactive = NULL;
112083c8fe1dSLisandro Dalcin   bnk->Hpre_inactive = NULL;
1121eb910715SAlp Dener   PetscFunctionReturn(0);
1122eb910715SAlp Dener }
1123eb910715SAlp Dener 
1124eb910715SAlp Dener /*------------------------------------------------------------*/
1125df278d8fSAlp Dener 
1126e0ed867bSAlp Dener PetscErrorCode TaoDestroy_BNK(Tao tao)
1127eb910715SAlp Dener {
1128eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1129eb910715SAlp Dener   PetscErrorCode ierr;
1130eb910715SAlp Dener 
1131eb910715SAlp Dener   PetscFunctionBegin;
1132eb910715SAlp Dener   if (tao->setupcalled) {
1133eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1134eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1135eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
113609164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
113709164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
113809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
113909164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
114062675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
114162675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1142c4b75bccSAlp Dener   }
1143ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1144ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1145ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1146ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1147ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
1148c4b75bccSAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1149c4b75bccSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1150ca964c49SAlp Dener   ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1151eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1152eb910715SAlp Dener   PetscFunctionReturn(0);
1153eb910715SAlp Dener }
1154eb910715SAlp Dener 
1155eb910715SAlp Dener /*------------------------------------------------------------*/
1156df278d8fSAlp Dener 
1157e0ed867bSAlp Dener PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1158eb910715SAlp Dener {
1159eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1160eb910715SAlp Dener   PetscErrorCode ierr;
1161eb910715SAlp Dener 
1162eb910715SAlp Dener   PetscFunctionBegin;
11634f4fdda4SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
116483c8fe1dSLisandro Dalcin   ierr = PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL);CHKERRQ(ierr);
116583c8fe1dSLisandro Dalcin   ierr = PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL);CHKERRQ(ierr);
116683c8fe1dSLisandro Dalcin   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL);CHKERRQ(ierr);
1167748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1168748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1169748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1170748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1171748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1172748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1173748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1174748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1175748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1176748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1177748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1178748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1179748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1180748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1181748696b1SAlp 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);
1182748696b1SAlp 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);
1183748696b1SAlp 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);
1184748696b1SAlp 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);
1185748696b1SAlp 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);
1186748696b1SAlp 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);
1187748696b1SAlp 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);
1188748696b1SAlp 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);
1189748696b1SAlp 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);
1190748696b1SAlp 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);
1191748696b1SAlp 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);
1192748696b1SAlp 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);
1193748696b1SAlp 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);
1194748696b1SAlp 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);
1195748696b1SAlp 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);
1196748696b1SAlp 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);
1197748696b1SAlp 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);
1198748696b1SAlp 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);
1199748696b1SAlp 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);
1200748696b1SAlp 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);
1201748696b1SAlp 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);
1202748696b1SAlp 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);
1203748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1204748696b1SAlp 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);
1205748696b1SAlp 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);
1206748696b1SAlp 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);
1207748696b1SAlp 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);
1208748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1209748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1210748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1211748696b1SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1212748696b1SAlp 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);
1213748696b1SAlp 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);
1214c0f10754SAlp 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);
1215eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1216*8ebe3e4eSStefano Zampini 
1217*8ebe3e4eSStefano Zampini   ierr = TaoSetOptionsPrefix(bnk->bncg,((PetscObject)(tao))->prefix);CHKERRQ(ierr);
1218*8ebe3e4eSStefano Zampini   ierr = TaoAppendOptionsPrefix(bnk->bncg,"tao_bnk_");CHKERRQ(ierr);
121933c78596SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1220*8ebe3e4eSStefano Zampini 
1221*8ebe3e4eSStefano Zampini   ierr = KSPSetOptionsPrefix(tao->ksp,((PetscObject)(tao))->prefix);CHKERRQ(ierr);
1222*8ebe3e4eSStefano Zampini   ierr = KSPAppendOptionsPrefix(tao->ksp,"tao_bnk_");CHKERRQ(ierr);
1223eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1224eb910715SAlp Dener   PetscFunctionReturn(0);
1225eb910715SAlp Dener }
1226eb910715SAlp Dener 
1227eb910715SAlp Dener /*------------------------------------------------------------*/
1228df278d8fSAlp Dener 
1229e0ed867bSAlp Dener PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1230eb910715SAlp Dener {
1231eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1232eb910715SAlp Dener   PetscInt       nrejects;
1233eb910715SAlp Dener   PetscBool      isascii;
1234eb910715SAlp Dener   PetscErrorCode ierr;
1235eb910715SAlp Dener 
1236eb910715SAlp Dener   PetscFunctionBegin;
1237eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1238eb910715SAlp Dener   if (isascii) {
1239eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1240b9ac7092SAlp Dener     if (bnk->M) {
1241cd929ea3SAlp Dener       ierr = MatLMVMGetRejectCount(bnk->M,&nrejects);CHKERRQ(ierr);
1242b9ac7092SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %D\n",nrejects);CHKERRQ(ierr);
1243eb910715SAlp Dener     }
1244e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1245eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1246b9ac7092SAlp Dener     if (bnk->M) {
1247eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1248b9ac7092SAlp Dener     }
1249eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1250eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1251eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1252eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1253eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1254eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1255eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1256eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1257eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1258eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1259eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1260eb910715SAlp Dener   }
1261eb910715SAlp Dener   PetscFunctionReturn(0);
1262eb910715SAlp Dener }
1263eb910715SAlp Dener 
1264eb910715SAlp Dener /* ---------------------------------------------------------- */
1265df278d8fSAlp Dener 
1266eb910715SAlp Dener /*MC
1267eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
126866ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1269eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1270eb910715SAlp Dener               Hk dk = -gk
12712b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12722b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1273eb910715SAlp Dener 
1274eb910715SAlp Dener     Options Database Keys:
12759fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
12769fa2b5dcSStefano Zampini . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
12779fa2b5dcSStefano Zampini . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
12789fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
12799fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
12809fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
12819fa2b5dcSStefano Zampini . -tao_bnk_sval - (developer) Hessian perturbation starting value
12829fa2b5dcSStefano Zampini . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
12839fa2b5dcSStefano Zampini . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
12849fa2b5dcSStefano Zampini . -tao_bnk_pmin - (developer) minimum Hessian perturbation
12859fa2b5dcSStefano Zampini . -tao_bnk_pmax - (developer) aximum Hessian perturbation
12869fa2b5dcSStefano Zampini . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
12879fa2b5dcSStefano Zampini . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
12889fa2b5dcSStefano Zampini . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
12899fa2b5dcSStefano Zampini . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
12909fa2b5dcSStefano Zampini . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
12919fa2b5dcSStefano Zampini . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction)
12929fa2b5dcSStefano Zampini . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
12939fa2b5dcSStefano Zampini . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
12949fa2b5dcSStefano Zampini . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction)
12959fa2b5dcSStefano Zampini . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
12969fa2b5dcSStefano Zampini . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
12979fa2b5dcSStefano Zampini . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
12989fa2b5dcSStefano Zampini . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
12999fa2b5dcSStefano Zampini . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
13009fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
13019fa2b5dcSStefano Zampini . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
13029fa2b5dcSStefano Zampini . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation)
13039fa2b5dcSStefano Zampini . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
13049fa2b5dcSStefano Zampini . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
13059fa2b5dcSStefano Zampini . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
13069fa2b5dcSStefano Zampini . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
13079fa2b5dcSStefano Zampini . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation)
13089fa2b5dcSStefano Zampini . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step)
13099fa2b5dcSStefano Zampini . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
13109fa2b5dcSStefano Zampini . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step)
13119fa2b5dcSStefano Zampini . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step)
13129fa2b5dcSStefano Zampini . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
13139fa2b5dcSStefano Zampini . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
13149fa2b5dcSStefano Zampini . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step)
13159fa2b5dcSStefano Zampini . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
13169fa2b5dcSStefano Zampini . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
13179fa2b5dcSStefano Zampini . -tao_bnk_mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
13189fa2b5dcSStefano Zampini . -tao_bnk_mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
13199fa2b5dcSStefano Zampini . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
13209fa2b5dcSStefano Zampini . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
13219fa2b5dcSStefano Zampini . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
13229fa2b5dcSStefano Zampini . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
13239fa2b5dcSStefano Zampini - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1324eb910715SAlp Dener 
1325eb910715SAlp Dener   Level: beginner
1326eb910715SAlp Dener M*/
1327eb910715SAlp Dener 
1328eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1329eb910715SAlp Dener {
1330eb910715SAlp Dener   TAO_BNK        *bnk;
1331eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1332eb910715SAlp Dener   PetscErrorCode ierr;
1333b9ac7092SAlp Dener   PC             pc;
1334eb910715SAlp Dener 
1335eb910715SAlp Dener   PetscFunctionBegin;
1336eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1337eb910715SAlp Dener 
1338eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1339eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1340eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1341eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1342eb910715SAlp Dener 
1343eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1344eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1345eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1346eb910715SAlp Dener 
1347eb910715SAlp Dener   tao->data = (void*)bnk;
1348eb910715SAlp Dener 
134966ed3702SAlp Dener   /*  Hessian shifting parameters */
1350e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1351e0ed867bSAlp Dener   bnk->computestep = TaoBNKComputeStep;
1352e0ed867bSAlp Dener 
1353eb910715SAlp Dener   bnk->sval   = 0.0;
1354eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1355eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1356eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1357eb910715SAlp Dener 
1358eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1359eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1360eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1361eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1362eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1363eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1364eb910715SAlp Dener 
1365eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1366eb910715SAlp Dener   bnk->nu1 = 0.25;
1367eb910715SAlp Dener   bnk->nu2 = 0.50;
1368eb910715SAlp Dener   bnk->nu3 = 1.00;
1369eb910715SAlp Dener   bnk->nu4 = 1.25;
1370eb910715SAlp Dener 
1371eb910715SAlp Dener   bnk->omega1 = 0.25;
1372eb910715SAlp Dener   bnk->omega2 = 0.50;
1373eb910715SAlp Dener   bnk->omega3 = 1.00;
1374eb910715SAlp Dener   bnk->omega4 = 2.00;
1375eb910715SAlp Dener   bnk->omega5 = 4.00;
1376eb910715SAlp Dener 
1377eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1378eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1379eb910715SAlp Dener   bnk->eta2 = 0.25;
1380eb910715SAlp Dener   bnk->eta3 = 0.50;
1381eb910715SAlp Dener   bnk->eta4 = 0.90;
1382eb910715SAlp Dener 
1383eb910715SAlp Dener   bnk->alpha1 = 0.25;
1384eb910715SAlp Dener   bnk->alpha2 = 0.50;
1385eb910715SAlp Dener   bnk->alpha3 = 1.00;
1386eb910715SAlp Dener   bnk->alpha4 = 2.00;
1387eb910715SAlp Dener   bnk->alpha5 = 4.00;
1388eb910715SAlp Dener 
1389eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1390eb910715SAlp Dener   bnk->mu1 = 0.10;
1391eb910715SAlp Dener   bnk->mu2 = 0.50;
1392eb910715SAlp Dener 
1393eb910715SAlp Dener   bnk->gamma1 = 0.25;
1394eb910715SAlp Dener   bnk->gamma2 = 0.50;
1395eb910715SAlp Dener   bnk->gamma3 = 2.00;
1396eb910715SAlp Dener   bnk->gamma4 = 4.00;
1397eb910715SAlp Dener 
1398eb910715SAlp Dener   bnk->theta = 0.05;
1399eb910715SAlp Dener 
1400eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1401eb910715SAlp Dener   bnk->mu1_i = 0.35;
1402eb910715SAlp Dener   bnk->mu2_i = 0.50;
1403eb910715SAlp Dener 
1404eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1405eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1406eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1407eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1408eb910715SAlp Dener 
1409eb910715SAlp Dener   bnk->theta_i = 0.25;
1410eb910715SAlp Dener 
1411eb910715SAlp Dener   /*  Remaining parameters */
1412c0f10754SAlp Dener   bnk->max_cg_its = 0;
1413eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1414eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1415770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14160a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14170a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
141862675beeSAlp Dener   bnk->dmin = 1.0e-6;
141962675beeSAlp Dener   bnk->dmax = 1.0e6;
1420eb910715SAlp Dener 
142183c8fe1dSLisandro Dalcin   bnk->M           = NULL;
142283c8fe1dSLisandro Dalcin   bnk->bfgs_pre    = NULL;
1423eb910715SAlp Dener   bnk->init_type   = BNK_INIT_INTERPOLATION;
14247b1c7716SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION;
14252f75a4aaSAlp Dener   bnk->as_type     = BNK_AS_BERTSEKAS;
1426eb910715SAlp Dener 
1427e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1428e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1429e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1430e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1431e031d6f5SAlp Dener 
1432c0f10754SAlp Dener   /* Create the line search */
1433eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1434eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1435eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1436eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1437eb910715SAlp Dener 
1438eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1439eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1440eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
144105de396fSBarry Smith   ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr);
1442f5a7d1c1SBarry Smith   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
1443b9ac7092SAlp Dener   ierr = PCSetType(pc, PCLMVM);CHKERRQ(ierr);
1444eb910715SAlp Dener   PetscFunctionReturn(0);
1445eb910715SAlp Dener }
1446