xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision f4db9bf761fa2b8a53a11d1c2fa5aeb4622be051)
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   TAO_BNK           *bnk = (TAO_BNK *)tao->data;
16eb910715SAlp Dener   PC                pc;
1789da521bSAlp Dener   PetscReal         f_min, ftrial, prered, actred, kappa, sigma, resnorm;
18eb910715SAlp Dener   PetscReal         tau, tau_1, tau_2, tau_max, tau_min, max_radius;
190ad3a497SAlp Dener   PetscBool         is_bfgs, is_jacobi, is_symmetric, sym_set;
20c4b75bccSAlp Dener   PetscInt          n, N, nDiff;
21eb910715SAlp Dener   PetscInt          i_max = 5;
22eb910715SAlp Dener   PetscInt          j_max = 1;
23eb910715SAlp Dener   PetscInt          i, j;
242e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
25eb910715SAlp Dener 
26eb910715SAlp Dener   PetscFunctionBegin;
2728017e9fSAlp Dener   /* Project the current point onto the feasible set */
289566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
299566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU));
30b9ac7092SAlp Dener   if (tao->bounded) {
319566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU));
32cd929ea3SAlp Dener   }
3328017e9fSAlp Dener 
3428017e9fSAlp Dener   /* Project the initial point onto the feasible region */
359566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution));
3628017e9fSAlp Dener 
3728017e9fSAlp Dener   /* Check convergence criteria */
389566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
399566063dSJacob Faibussowitsch   PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
409566063dSJacob Faibussowitsch   PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
419566063dSJacob Faibussowitsch   PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
429566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
4328017e9fSAlp Dener 
44c0f10754SAlp Dener   /* Test the initial point for convergence */
459566063dSJacob Faibussowitsch   PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
469566063dSJacob Faibussowitsch   PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
473c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
489566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its));
499566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0));
509566063dSJacob Faibussowitsch   PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
51c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
52c0f10754SAlp Dener 
53e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
54eb910715SAlp Dener   bnk->ksp_atol = 0;
55eb910715SAlp Dener   bnk->ksp_rtol = 0;
56eb910715SAlp Dener   bnk->ksp_dtol = 0;
57eb910715SAlp Dener   bnk->ksp_ctol = 0;
58eb910715SAlp Dener   bnk->ksp_negc = 0;
59eb910715SAlp Dener   bnk->ksp_iter = 0;
60eb910715SAlp Dener   bnk->ksp_othr = 0;
61eb910715SAlp Dener 
62e031d6f5SAlp Dener   /* Reset accepted step type counters */
63e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
64e031d6f5SAlp Dener   bnk->newt = 0;
65e031d6f5SAlp Dener   bnk->bfgs = 0;
66e031d6f5SAlp Dener   bnk->sgrad = 0;
67e031d6f5SAlp Dener   bnk->grad = 0;
68e031d6f5SAlp Dener 
69fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
70fed79b8eSAlp Dener   bnk->pert = bnk->sval;
71fed79b8eSAlp Dener 
72937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
739566063dSJacob Faibussowitsch   PetscCall(VecSet(tao->stepdirection, 0.0));
74937a31a1SAlp Dener 
75e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
769566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
79b9ac7092SAlp Dener   if (is_bfgs) {
80b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
819566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M));
829566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
839566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
849566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(bnk->M, n, n, N, N));
859566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient));
869566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric));
873c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
88b9ac7092SAlp Dener   } else if (is_jacobi) {
899566063dSJacob Faibussowitsch     PetscCall(PCJacobiSetUseAbs(pc,PETSC_TRUE));
90e031d6f5SAlp Dener   }
91e031d6f5SAlp Dener 
92e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
939566063dSJacob Faibussowitsch   PetscCall(VecSet(bnk->Diag_min, bnk->dmin));
949566063dSJacob Faibussowitsch   PetscCall(VecSet(bnk->Diag_max, bnk->dmax));
95eb910715SAlp Dener 
96eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
97eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
98c0f10754SAlp Dener   *needH = PETSC_TRUE;
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp,"KSPCGSetRadius_C",&kspTR));
1002e6e4ca1SStefano Zampini   if (kspTR) {
10162675beeSAlp Dener     switch (initType) {
102eb910715SAlp Dener     case BNK_INIT_CONSTANT:
103eb910715SAlp Dener       /* Use the initial radius specified */
104c0f10754SAlp Dener       tao->trust = tao->trust0;
105eb910715SAlp Dener       break;
106eb910715SAlp Dener 
107eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
108c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
109eb910715SAlp Dener       max_radius = 0.0;
11008752603SAlp Dener       tao->trust = tao->trust0;
111eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1120a4511e9SAlp Dener         f_min = bnk->f;
113eb910715SAlp Dener         sigma = 0.0;
114eb910715SAlp Dener 
115c0f10754SAlp Dener         if (*needH) {
11662602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
1179566063dSJacob Faibussowitsch           PetscCall((*bnk->computehessian)(tao));
1189566063dSJacob Faibussowitsch           PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE));
1199566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&bnk->H_inactive));
12089da521bSAlp Dener           if (bnk->active_idx) {
1219566063dSJacob Faibussowitsch             PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
12228017e9fSAlp Dener           } else {
1239566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)tao->hessian));
124c5e9d94cSAlp Dener             bnk->H_inactive = tao->hessian;
12528017e9fSAlp Dener           }
126c0f10754SAlp Dener           *needH = PETSC_FALSE;
127eb910715SAlp Dener         }
128eb910715SAlp Dener 
129eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
13062602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
1319566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->Xold));
1329566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient));
1339566063dSJacob Faibussowitsch           PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution));
13489da521bSAlp Dener           /* Compute the step we actually accepted */
1359566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->W));
1369566063dSJacob Faibussowitsch           PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold));
13762602cfbSAlp Dener           /* Compute the objective at the trial */
1389566063dSJacob Faibussowitsch           PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial));
1393c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(bnk->f),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
1409566063dSJacob Faibussowitsch           PetscCall(VecCopy(bnk->Xold, tao->solution));
141eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
142eb910715SAlp Dener             tau = bnk->gamma1_i;
143eb910715SAlp Dener           } else {
1440a4511e9SAlp Dener             if (ftrial < f_min) {
1450a4511e9SAlp Dener               f_min = ftrial;
146eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
147eb910715SAlp Dener             }
14808752603SAlp Dener 
149770b7498SAlp Dener             /* Compute the predicted and actual reduction */
15089da521bSAlp Dener             if (bnk->active_idx) {
1519566063dSJacob Faibussowitsch               PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1529566063dSJacob Faibussowitsch               PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1532ab2a32cSAlp Dener             } else {
15408752603SAlp Dener               bnk->X_inactive = bnk->W;
15508752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1562ab2a32cSAlp Dener             }
1579566063dSJacob Faibussowitsch             PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
1589566063dSJacob Faibussowitsch             PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered));
15989da521bSAlp Dener             if (bnk->active_idx) {
1609566063dSJacob Faibussowitsch               PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1619566063dSJacob Faibussowitsch               PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1622ab2a32cSAlp Dener             }
163eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
164eb910715SAlp Dener             actred = bnk->f - ftrial;
1653105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
166eb910715SAlp Dener               kappa = 1.0;
1673105154fSTodd Munson             } else {
168eb910715SAlp Dener               kappa = actred / prered;
169eb910715SAlp Dener             }
170eb910715SAlp Dener 
171eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
172eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
173eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
174eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
175eb910715SAlp Dener 
17618cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
177eb910715SAlp Dener               /*  Great agreement */
178eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
179eb910715SAlp Dener 
180eb910715SAlp Dener               if (tau_max < 1.0) {
181eb910715SAlp Dener                 tau = bnk->gamma3_i;
1823105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
183eb910715SAlp Dener                 tau = bnk->gamma4_i;
1843105154fSTodd Munson               } else {
185eb910715SAlp Dener                 tau = tau_max;
186eb910715SAlp Dener               }
18718cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
188eb910715SAlp Dener               /*  Good agreement */
189eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
190eb910715SAlp Dener 
191eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
192eb910715SAlp Dener                 tau = bnk->gamma2_i;
193eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
194eb910715SAlp Dener                 tau = bnk->gamma3_i;
195eb910715SAlp Dener               } else {
196eb910715SAlp Dener                 tau = tau_max;
197eb910715SAlp Dener               }
1988f8a4e06SAlp Dener             } else {
199eb910715SAlp Dener               /*  Not good agreement */
200eb910715SAlp Dener               if (tau_min > 1.0) {
201eb910715SAlp Dener                 tau = bnk->gamma2_i;
202eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
203eb910715SAlp Dener                 tau = bnk->gamma1_i;
204eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
205eb910715SAlp Dener                 tau = bnk->gamma1_i;
2063105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
207eb910715SAlp Dener                 tau = tau_1;
2083105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
209eb910715SAlp Dener                 tau = tau_2;
210eb910715SAlp Dener               } else {
211eb910715SAlp Dener                 tau = tau_max;
212eb910715SAlp Dener               }
213eb910715SAlp Dener             }
214eb910715SAlp Dener           }
215eb910715SAlp Dener           tao->trust = tau * tao->trust;
216eb910715SAlp Dener         }
217eb910715SAlp Dener 
2180a4511e9SAlp Dener         if (f_min < bnk->f) {
219937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2200a4511e9SAlp Dener           bnk->f = f_min;
2219566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->Xold));
2229566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution,sigma,tao->gradient));
2239566063dSJacob Faibussowitsch           PetscCall(TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution));
2249566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, tao->stepdirection));
2259566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
2269566063dSJacob Faibussowitsch           PetscCall(TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient));
2279566063dSJacob Faibussowitsch           PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
2289566063dSJacob Faibussowitsch           PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
2299566063dSJacob Faibussowitsch           PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
230937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
2319566063dSJacob Faibussowitsch           PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
232c0f10754SAlp Dener           *needH = PETSC_TRUE;
233937a31a1SAlp Dener           /* Test the new step for convergence */
2349566063dSJacob Faibussowitsch           PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
2359566063dSJacob Faibussowitsch           PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
2363c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(resnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2379566063dSJacob Faibussowitsch           PetscCall(TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its));
2389566063dSJacob Faibussowitsch           PetscCall(TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0));
2399566063dSJacob Faibussowitsch           PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP));
240eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
241937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
2429566063dSJacob Faibussowitsch           PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
243eb910715SAlp Dener         }
244eb910715SAlp Dener       }
245eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
246e031d6f5SAlp Dener 
247e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
248e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
249e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
250eb910715SAlp Dener       break;
251eb910715SAlp Dener 
252eb910715SAlp Dener     default:
253eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
254eb910715SAlp Dener       tao->trust = 0.0;
255eb910715SAlp Dener       break;
256eb910715SAlp Dener     }
257eb910715SAlp Dener   }
258eb910715SAlp Dener   PetscFunctionReturn(0);
259eb910715SAlp Dener }
260eb910715SAlp Dener 
261df278d8fSAlp Dener /*------------------------------------------------------------*/
262df278d8fSAlp Dener 
263e0ed867bSAlp Dener /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */
26462675beeSAlp Dener 
26562675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
26662675beeSAlp Dener {
26762675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
26862675beeSAlp Dener 
26962675beeSAlp Dener   PetscFunctionBegin;
27062675beeSAlp Dener   /* Compute the Hessian */
2719566063dSJacob Faibussowitsch   PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre));
27262675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
273b9ac7092SAlp Dener   if (bnk->M) {
2749566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
27562675beeSAlp Dener   }
276e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
2779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
2789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
279f5766c09SAlp Dener   if (bnk->active_idx) {
2809566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
281e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
2829566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
283e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
284e0ed867bSAlp Dener     } else {
2859566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive));
286e0ed867bSAlp Dener     }
287e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
2889566063dSJacob Faibussowitsch       PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx));
289e0ed867bSAlp Dener     }
290e0ed867bSAlp Dener   } else {
2919566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)tao->hessian));
292c5e9d94cSAlp Dener     bnk->H_inactive = tao->hessian;
293e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
2949566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
295e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
296e0ed867bSAlp Dener     } else {
2979566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre));
298c5e9d94cSAlp Dener       bnk->Hpre_inactive = tao->hessian_pre;
299e0ed867bSAlp Dener     }
300e0ed867bSAlp Dener     if (bnk->bfgs_pre) {
3019566063dSJacob Faibussowitsch       PetscCall(PCLMVMClearIS(bnk->bfgs_pre));
302e0ed867bSAlp Dener     }
303e0ed867bSAlp Dener   }
30462675beeSAlp Dener   PetscFunctionReturn(0);
30562675beeSAlp Dener }
30662675beeSAlp Dener 
30762675beeSAlp Dener /*------------------------------------------------------------*/
30862675beeSAlp Dener 
3092f75a4aaSAlp Dener /* Routine for estimating the active set */
3102f75a4aaSAlp Dener 
31108752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3122f75a4aaSAlp Dener {
3132f75a4aaSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
314*f4db9bf7SStefano Zampini   PetscBool      hessComputed, diagExists, hadactive;
3152f75a4aaSAlp Dener 
3162f75a4aaSAlp Dener   PetscFunctionBegin;
317*f4db9bf7SStefano Zampini   hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE;
31808752603SAlp Dener   switch (asType) {
3192f75a4aaSAlp Dener   case BNK_AS_NONE:
3209566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bnk->inactive_idx));
3219566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx));
3229566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bnk->active_idx));
3239566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx));
3242f75a4aaSAlp Dener     break;
3252f75a4aaSAlp Dener 
3262f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3272f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
328b9ac7092SAlp Dener     if (bnk->M) {
3292f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3309566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W));
3312f75a4aaSAlp Dener     } else {
332fc5ca067SStefano Zampini       hessComputed = diagExists = PETSC_FALSE;
333f5766c09SAlp Dener       if (tao->hessian) {
3349566063dSJacob Faibussowitsch         PetscCall(MatAssembled(tao->hessian, &hessComputed));
335f5766c09SAlp Dener       }
336fc5ca067SStefano Zampini       if (hessComputed) {
3379566063dSJacob Faibussowitsch         PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists));
338fc5ca067SStefano Zampini       }
339fc5ca067SStefano Zampini       if (diagExists) {
3409b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3419566063dSJacob Faibussowitsch         PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork));
3429566063dSJacob Faibussowitsch         PetscCall(VecAbs(bnk->Xwork));
3439566063dSJacob Faibussowitsch         PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork));
3449566063dSJacob Faibussowitsch         PetscCall(VecReciprocal(bnk->Xwork));
3459566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient));
34661be54a6SAlp Dener       } else {
347c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
3489566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W));
34961be54a6SAlp Dener       }
3502f75a4aaSAlp Dener     }
3519566063dSJacob Faibussowitsch     PetscCall(VecScale(bnk->W, -1.0));
352d0609cedSBarry Smith     PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
353d0609cedSBarry Smith                                       &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx));
354c4b75bccSAlp Dener     break;
3552f75a4aaSAlp Dener 
3562f75a4aaSAlp Dener   default:
3572f75a4aaSAlp Dener     break;
3582f75a4aaSAlp Dener   }
359*f4db9bf7SStefano Zampini   bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */
3602f75a4aaSAlp Dener   PetscFunctionReturn(0);
3612f75a4aaSAlp Dener }
3622f75a4aaSAlp Dener 
3632f75a4aaSAlp Dener /*------------------------------------------------------------*/
3642f75a4aaSAlp Dener 
3652f75a4aaSAlp Dener /* Routine for bounding the step direction */
3662f75a4aaSAlp Dener 
367a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
3682f75a4aaSAlp Dener {
3692f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3702f75a4aaSAlp Dener 
3712f75a4aaSAlp Dener   PetscFunctionBegin;
372a1318120SAlp Dener   switch (asType) {
3732f75a4aaSAlp Dener   case BNK_AS_NONE:
3749566063dSJacob Faibussowitsch     PetscCall(VecISSet(step, bnk->active_idx, 0.0));
3752f75a4aaSAlp Dener     break;
3762f75a4aaSAlp Dener 
3772f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3789566063dSJacob Faibussowitsch     PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step));
3792f75a4aaSAlp Dener     break;
3802f75a4aaSAlp Dener 
3812f75a4aaSAlp Dener   default:
3822f75a4aaSAlp Dener     break;
3832f75a4aaSAlp Dener   }
3842f75a4aaSAlp Dener   PetscFunctionReturn(0);
3852f75a4aaSAlp Dener }
3862f75a4aaSAlp Dener 
387e031d6f5SAlp Dener /*------------------------------------------------------------*/
388e031d6f5SAlp Dener 
389e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
390e031d6f5SAlp Dener    accelerate Newton convergence.
391e031d6f5SAlp Dener 
392e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
393e031d6f5SAlp Dener    for more gradient evaluations.
394e031d6f5SAlp Dener */
395e031d6f5SAlp Dener 
396c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
397c0f10754SAlp Dener {
398c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
399c0f10754SAlp Dener 
400c0f10754SAlp Dener   PetscFunctionBegin;
401c0f10754SAlp Dener   *terminate = PETSC_FALSE;
402c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
403c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
404c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
405c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
4069566063dSJacob Faibussowitsch     PetscCall(TaoSolve(bnk->bncg));
407c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
408c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
409c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
410c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
411c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
412e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
413c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
414c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
415c0f10754SAlp 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) {
416c0f10754SAlp Dener       *terminate = PETSC_TRUE;
41761be54a6SAlp Dener     } else {
4189566063dSJacob Faibussowitsch       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
419c0f10754SAlp Dener     }
420c0f10754SAlp Dener   }
421c0f10754SAlp Dener   PetscFunctionReturn(0);
422c0f10754SAlp Dener }
423c0f10754SAlp Dener 
4242f75a4aaSAlp Dener /*------------------------------------------------------------*/
4252f75a4aaSAlp Dener 
426c0f10754SAlp Dener /* Routine for computing the Newton step. */
427df278d8fSAlp Dener 
4286b591159SAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
429eb910715SAlp Dener {
430eb910715SAlp Dener   TAO_BNK           *bnk = (TAO_BNK *)tao->data;
431eb910715SAlp Dener   PetscInt          bfgsUpdates = 0;
432eb910715SAlp Dener   PetscInt          kspits;
433bddd1ffdSAlp Dener   PetscBool         is_lmvm;
4342e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
435eb910715SAlp Dener 
436eb910715SAlp Dener   PetscFunctionBegin;
43789da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
43889da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
43989da521bSAlp Dener   if (!bnk->inactive_idx) {
4409566063dSJacob Faibussowitsch     PetscCall(VecSet(tao->stepdirection, 0.0));
4419566063dSJacob Faibussowitsch     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
44289da521bSAlp Dener     PetscFunctionReturn(0);
44389da521bSAlp Dener   }
44489da521bSAlp Dener 
44562675beeSAlp Dener   /* Shift the reduced Hessian matrix */
446e831869dSStefano Zampini   if (shift && bnk->pert > 0) {
4479566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm));
448f7bf01afSAlp Dener     if (is_lmvm) {
4499566063dSJacob Faibussowitsch       PetscCall(MatShift(tao->hessian, bnk->pert));
450f7bf01afSAlp Dener     } else {
4519566063dSJacob Faibussowitsch       PetscCall(MatShift(bnk->H_inactive, bnk->pert));
45262675beeSAlp Dener       if (bnk->H_inactive != bnk->Hpre_inactive) {
4539566063dSJacob Faibussowitsch         PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert));
45462675beeSAlp Dener       }
45562675beeSAlp Dener     }
456f7bf01afSAlp Dener   }
45762675beeSAlp Dener 
458eb910715SAlp Dener   /* Solve the Newton system of equations */
459937a31a1SAlp Dener   tao->ksp_its = 0;
4609566063dSJacob Faibussowitsch   PetscCall(VecSet(tao->stepdirection, 0.0));
461*f4db9bf7SStefano Zampini   if (bnk->resetksp) {
4629566063dSJacob Faibussowitsch     PetscCall(KSPReset(tao->ksp));
4639566063dSJacob Faibussowitsch     PetscCall(KSPResetFromOptions(tao->ksp));
464*f4db9bf7SStefano Zampini     bnk->resetksp = PETSC_FALSE;
465*f4db9bf7SStefano Zampini   }
4669566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive));
4679566063dSJacob Faibussowitsch   PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork));
46889da521bSAlp Dener   if (bnk->active_idx) {
4699566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4709566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4715e9b73cbSAlp Dener   } else {
4725e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4735e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
47428017e9fSAlp Dener   }
4759566063dSJacob Faibussowitsch   PetscCall(KSPCGSetRadius(tao->ksp,tao->trust));
4769566063dSJacob Faibussowitsch   PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4779566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(tao->ksp,&kspits));
478eb910715SAlp Dener   tao->ksp_its += kspits;
479eb910715SAlp Dener   tao->ksp_tot_its += kspits;
480*f4db9bf7SStefano Zampini   PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp,"KSPCGGetNormD_C",&kspTR));
481*f4db9bf7SStefano Zampini   if (kspTR) {
4829566063dSJacob Faibussowitsch     PetscCall(KSPCGGetNormD(tao->ksp,&bnk->dnorm));
483eb910715SAlp Dener 
484eb910715SAlp Dener     if (0.0 == tao->trust) {
485eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
486080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
487080d2917SAlp Dener         tao->trust = bnk->dnorm;
488eb910715SAlp Dener 
489eb910715SAlp Dener         /* Modify the radius if it is too large or small */
490eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
491eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
492eb910715SAlp Dener       } else {
493eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
494eb910715SAlp Dener            the trust-region subproblem to get a direction */
495eb910715SAlp Dener         tao->trust = tao->trust0;
496eb910715SAlp Dener 
497eb910715SAlp Dener         /* Modify the radius if it is too large or small */
498eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
499eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
500eb910715SAlp Dener 
5019566063dSJacob Faibussowitsch         PetscCall(KSPCGSetRadius(tao->ksp,tao->trust));
5029566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
5039566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp,&kspits));
504eb910715SAlp Dener         tao->ksp_its += kspits;
505eb910715SAlp Dener         tao->ksp_tot_its += kspits;
5069566063dSJacob Faibussowitsch         PetscCall(KSPCGGetNormD(tao->ksp,&bnk->dnorm));
507eb910715SAlp Dener 
5083c859ba3SBarry Smith         PetscCheck(bnk->dnorm != 0.0,PetscObjectComm((PetscObject)tao),PETSC_ERR_PLIB, "Initial direction zero");
509eb910715SAlp Dener       }
510eb910715SAlp Dener     }
511eb910715SAlp Dener   }
5125e9b73cbSAlp Dener   /* Restore sub vectors back */
51389da521bSAlp Dener   if (bnk->active_idx) {
5149566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5159566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5165e9b73cbSAlp Dener   }
517770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
5189566063dSJacob Faibussowitsch   PetscCall(VecScale(tao->stepdirection, -1.0));
5199566063dSJacob Faibussowitsch   PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
520770b7498SAlp Dener 
521770b7498SAlp Dener   /* Record convergence reasons */
5229566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason));
523e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
524770b7498SAlp Dener     ++bnk->ksp_atol;
525e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
526770b7498SAlp Dener     ++bnk->ksp_rtol;
527e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
528770b7498SAlp Dener     ++bnk->ksp_ctol;
529e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
530770b7498SAlp Dener     ++bnk->ksp_negc;
531e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
532770b7498SAlp Dener     ++bnk->ksp_dtol;
533e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
534770b7498SAlp Dener     ++bnk->ksp_iter;
535770b7498SAlp Dener   } else {
536770b7498SAlp Dener     ++bnk->ksp_othr;
537770b7498SAlp Dener   }
538fed79b8eSAlp Dener 
539fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
540b9ac7092SAlp Dener   if (bnk->M) {
5419566063dSJacob Faibussowitsch     PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
542b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
543fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
5449566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
5459566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
546eb910715SAlp Dener     }
547fed79b8eSAlp Dener   }
5486b591159SAlp Dener   *step_type = BNK_NEWTON;
549e465cd6fSAlp Dener   PetscFunctionReturn(0);
550e465cd6fSAlp Dener }
551eb910715SAlp Dener 
55262675beeSAlp Dener /*------------------------------------------------------------*/
55362675beeSAlp Dener 
5545e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5555e9b73cbSAlp Dener 
5565e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
5575e9b73cbSAlp Dener {
5585e9b73cbSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
5595e9b73cbSAlp Dener 
5605e9b73cbSAlp Dener   PetscFunctionBegin;
5615e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
56289da521bSAlp Dener   if (bnk->active_idx) {
5639566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5649566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5659566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5665e9b73cbSAlp Dener   } else {
5675e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
5685e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5695e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
5705e9b73cbSAlp Dener   }
5715e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5729566063dSJacob Faibussowitsch   PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
5739566063dSJacob Faibussowitsch   PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive));
5749566063dSJacob Faibussowitsch   PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered));
5755e9b73cbSAlp Dener   /* Restore the sub vectors */
57689da521bSAlp Dener   if (bnk->active_idx) {
5779566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5789566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5799566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5805e9b73cbSAlp Dener   }
5815e9b73cbSAlp Dener   PetscFunctionReturn(0);
5825e9b73cbSAlp Dener }
5835e9b73cbSAlp Dener 
5845e9b73cbSAlp Dener /*------------------------------------------------------------*/
5855e9b73cbSAlp Dener 
58662675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
58762675beeSAlp Dener 
58862675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
58962675beeSAlp Dener    in the event that the Newton step fails the test.
59062675beeSAlp Dener */
59162675beeSAlp Dener 
592e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
593e465cd6fSAlp Dener {
594e465cd6fSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
595b2d8c577SAlp Dener   PetscReal      gdx, e_min;
596e465cd6fSAlp Dener   PetscInt       bfgsUpdates;
597e465cd6fSAlp Dener 
598e465cd6fSAlp Dener   PetscFunctionBegin;
5996b591159SAlp Dener   switch (*stepType) {
6006b591159SAlp Dener   case BNK_NEWTON:
6019566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
602eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
603eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
604eb910715SAlp Dener         Update the perturbation for next time */
605eb910715SAlp Dener       if (bnk->pert <= 0.0) {
6062e6e4ca1SStefano Zampini         PetscBool is_gltr;
6072e6e4ca1SStefano Zampini 
608eb910715SAlp Dener         /* Initialize the perturbation */
609eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
6109566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr));
6112e6e4ca1SStefano Zampini         if (is_gltr) {
6129566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min));
613eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
614eb910715SAlp Dener         }
615eb910715SAlp Dener       } else {
616eb910715SAlp Dener         /* Increase the perturbation */
617eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
618eb910715SAlp Dener       }
619eb910715SAlp Dener 
6200ad3a497SAlp Dener       if (!bnk->M) {
621eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
622eb910715SAlp Dener           Must use gradient direction in this case */
6239566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->gradient, tao->stepdirection));
624eb910715SAlp Dener         *stepType = BNK_GRADIENT;
625eb910715SAlp Dener       } else {
626eb910715SAlp Dener         /* Attempt to use the BFGS direction */
6279566063dSJacob Faibussowitsch         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
628eb910715SAlp Dener 
6298d5ead36SAlp Dener         /* Check for success (descent direction)
6308d5ead36SAlp Dener           NOTE: Negative gdx here means not a descent direction because
6318d5ead36SAlp Dener           the fall-back step is missing a negative sign. */
6329566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
6333105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
634eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
635eb910715SAlp Dener             We can assert bfgsUpdates > 1 in this case because
636eb910715SAlp Dener             the first solve produces the scaled gradient direction,
637eb910715SAlp Dener             which is guaranteed to be descent */
638eb910715SAlp Dener 
639eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
6409566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6419566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6429566063dSJacob Faibussowitsch           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
643eb910715SAlp Dener 
644eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
645eb910715SAlp Dener         } else {
6469566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
647eb910715SAlp Dener           if (1 == bfgsUpdates) {
648eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
649eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
650eb910715SAlp Dener           } else {
651eb910715SAlp Dener             *stepType = BNK_BFGS;
652eb910715SAlp Dener           }
653eb910715SAlp Dener         }
654eb910715SAlp Dener       }
6558d5ead36SAlp Dener       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6569566063dSJacob Faibussowitsch       PetscCall(VecScale(tao->stepdirection, -1.0));
6579566063dSJacob Faibussowitsch       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
658eb910715SAlp Dener     } else {
659eb910715SAlp Dener       /* Computed Newton step is descent */
660eb910715SAlp Dener       switch (ksp_reason) {
661eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
662eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
663eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
664eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
665eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
666eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
667eb910715SAlp Dener         if (bnk->pert <= 0.0) {
6682e6e4ca1SStefano Zampini           PetscBool is_gltr;
6692e6e4ca1SStefano Zampini 
670eb910715SAlp Dener           /* Initialize the perturbation */
671eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
6729566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr));
6732e6e4ca1SStefano Zampini           if (is_gltr) {
6749566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
675eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
676eb910715SAlp Dener           }
677eb910715SAlp Dener         } else {
678eb910715SAlp Dener           /* Increase the perturbation */
679eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
680eb910715SAlp Dener         }
681eb910715SAlp Dener         break;
682eb910715SAlp Dener 
683eb910715SAlp Dener       default:
684eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
685eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
686eb910715SAlp Dener         if (bnk->pert < bnk->pmin) {
687eb910715SAlp Dener           bnk->pert = 0.0;
688eb910715SAlp Dener         }
689eb910715SAlp Dener         break;
690eb910715SAlp Dener       }
691fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
692eb910715SAlp Dener     }
6936b591159SAlp Dener     break;
6946b591159SAlp Dener 
6956b591159SAlp Dener   case BNK_BFGS:
6966b591159SAlp Dener     /* Check for success (descent direction) */
6979566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
6986b591159SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
6996b591159SAlp Dener       /* Step is not descent or solve was not successful
7006b591159SAlp Dener          Use steepest descent direction (scaled) */
7019566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7029566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7039566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection));
7049566063dSJacob Faibussowitsch       PetscCall(VecScale(tao->stepdirection,-1.0));
7059566063dSJacob Faibussowitsch       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
7066b591159SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
7076b591159SAlp Dener     } else {
7086b591159SAlp Dener       *stepType = BNK_BFGS;
7096b591159SAlp Dener     }
7106b591159SAlp Dener     break;
7116b591159SAlp Dener 
7126b591159SAlp Dener   case BNK_SCALED_GRADIENT:
7136b591159SAlp Dener     break;
7146b591159SAlp Dener 
7156b591159SAlp Dener   default:
7166b591159SAlp Dener     break;
7176b591159SAlp Dener   }
7186b591159SAlp Dener 
719eb910715SAlp Dener   PetscFunctionReturn(0);
720eb910715SAlp Dener }
721eb910715SAlp Dener 
722df278d8fSAlp Dener /*------------------------------------------------------------*/
723df278d8fSAlp Dener 
724df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
725df278d8fSAlp Dener 
726df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
727df278d8fSAlp Dener   Newton step does not produce a valid step length.
728df278d8fSAlp Dener */
729df278d8fSAlp Dener 
730937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
731c14b763aSAlp Dener {
732c14b763aSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
733c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
734b2d8c577SAlp Dener   PetscReal                    e_min, gdx;
735c14b763aSAlp Dener   PetscInt                     bfgsUpdates;
736c14b763aSAlp Dener 
737c14b763aSAlp Dener   PetscFunctionBegin;
738c14b763aSAlp Dener   /* Perform the linesearch */
7399566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7409566063dSJacob Faibussowitsch   PetscCall(TaoAddLineSearchCounts(tao));
741c14b763aSAlp Dener 
742b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
743c14b763aSAlp Dener     /* Linesearch failed, revert solution */
744c14b763aSAlp Dener     bnk->f = bnk->fold;
7459566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->Xold, tao->solution));
7469566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
747c14b763aSAlp Dener 
748937a31a1SAlp Dener     switch(*stepType) {
749c14b763aSAlp Dener     case BNK_NEWTON:
7508d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
751c14b763aSAlp Dener          Update the perturbation for next time */
752c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
7532e6e4ca1SStefano Zampini         PetscBool is_gltr;
7542e6e4ca1SStefano Zampini 
755c14b763aSAlp Dener         /* Initialize the perturbation */
756c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
7579566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr));
7582e6e4ca1SStefano Zampini         if (is_gltr) {
7599566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp,&e_min));
760c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
761c14b763aSAlp Dener         }
762c14b763aSAlp Dener       } else {
763c14b763aSAlp Dener         /* Increase the perturbation */
764c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
765c14b763aSAlp Dener       }
766c14b763aSAlp Dener 
7670ad3a497SAlp Dener       if (!bnk->M) {
768c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
769c14b763aSAlp Dener            Must use gradient direction in this case */
7709566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection));
771937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
772c14b763aSAlp Dener       } else {
773c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
7749566063dSJacob Faibussowitsch         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
7758d5ead36SAlp Dener         /* Check for success (descent direction)
7768d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
7779566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
7783105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
779c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
780c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
781c14b763aSAlp Dener              Use steepest descent direction (scaled) */
7829566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7839566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7849566063dSJacob Faibussowitsch           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
785c14b763aSAlp Dener 
786c14b763aSAlp Dener           bfgsUpdates = 1;
787937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
788c14b763aSAlp Dener         } else {
7899566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
790c14b763aSAlp Dener           if (1 == bfgsUpdates) {
791c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
792937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
793c14b763aSAlp Dener           } else {
794937a31a1SAlp Dener             *stepType = BNK_BFGS;
795c14b763aSAlp Dener           }
796c14b763aSAlp Dener         }
797c14b763aSAlp Dener       }
798c14b763aSAlp Dener       break;
799c14b763aSAlp Dener 
800c14b763aSAlp Dener     case BNK_BFGS:
801c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
802c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
803c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
8049566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
8059566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
8069566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
807c14b763aSAlp Dener 
808c14b763aSAlp Dener       bfgsUpdates = 1;
809937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
810c14b763aSAlp Dener       break;
811c14b763aSAlp Dener     }
8128d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8139566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
8149566063dSJacob Faibussowitsch     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
815c14b763aSAlp Dener 
8168d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
8179566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
8189566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
819c14b763aSAlp Dener   }
820c14b763aSAlp Dener   *reason = ls_reason;
821c14b763aSAlp Dener   PetscFunctionReturn(0);
822c14b763aSAlp Dener }
823c14b763aSAlp Dener 
824df278d8fSAlp Dener /*------------------------------------------------------------*/
825df278d8fSAlp Dener 
826df278d8fSAlp Dener /* Routine for updating the trust radius.
827df278d8fSAlp Dener 
828df278d8fSAlp Dener   Function features three different update methods:
829df278d8fSAlp Dener   1) Line-search step length based
830df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
831df278d8fSAlp Dener   3) Interpolation
832df278d8fSAlp Dener */
833df278d8fSAlp Dener 
83428017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
835080d2917SAlp Dener {
836080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
837080d2917SAlp Dener 
838b1c2d0e3SAlp Dener   PetscReal      step, kappa;
839080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
840080d2917SAlp Dener 
841080d2917SAlp Dener   PetscFunctionBegin;
842080d2917SAlp Dener   /* Update trust region radius */
843080d2917SAlp Dener   *accept = PETSC_FALSE;
84428017e9fSAlp Dener   switch(updateType) {
845080d2917SAlp Dener   case BNK_UPDATE_STEP:
846c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
847080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
8489566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step));
849080d2917SAlp Dener       if (step < bnk->nu1) {
850080d2917SAlp Dener         /* Very bad step taken; reduce radius */
851080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
852080d2917SAlp Dener       } else if (step < bnk->nu2) {
853080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
854080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
855080d2917SAlp Dener       } else if (step < bnk->nu3) {
856080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
857080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
858080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
859080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
860080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
861080d2917SAlp Dener         }
862080d2917SAlp Dener       } else if (step < bnk->nu4) {
863080d2917SAlp Dener         /*  Full step taken; increase the radius */
864080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
865080d2917SAlp Dener       } else {
866080d2917SAlp Dener         /*  More than full step taken; increase the radius */
867080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
868080d2917SAlp Dener       }
869080d2917SAlp Dener     } else {
870080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
871080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
872080d2917SAlp Dener     }
873080d2917SAlp Dener     break;
874080d2917SAlp Dener 
875080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
876080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
877e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
878fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
879fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
880fed79b8eSAlp Dener            be rejected! */
881080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8823105154fSTodd Munson       } else {
883b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
884080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
885080d2917SAlp Dener         } else {
8863105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
887080d2917SAlp Dener             kappa = 1.0;
8883105154fSTodd Munson           } else {
889080d2917SAlp Dener             kappa = actred / prered;
890080d2917SAlp Dener           }
891fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
892080d2917SAlp Dener           if (kappa < bnk->eta1) {
893fed79b8eSAlp Dener             /* Reject the step */
894080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8953105154fSTodd Munson           } else {
896fed79b8eSAlp Dener             /* Accept the step */
897c133c014SAlp Dener             *accept = PETSC_TRUE;
898c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8998d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
900080d2917SAlp Dener               if (kappa < bnk->eta2) {
901080d2917SAlp Dener                 /* Marginal bad step */
902c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
9033105154fSTodd Munson               } else if (kappa < bnk->eta3) {
904fed79b8eSAlp Dener                 /* Reasonable step */
905fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
9063105154fSTodd Munson               } else if (kappa < bnk->eta4) {
907080d2917SAlp Dener                 /* Good step */
908c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
9093105154fSTodd Munson               } else {
910080d2917SAlp Dener                 /* Very good step */
911c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
912080d2917SAlp Dener               }
913c133c014SAlp Dener             }
914080d2917SAlp Dener           }
915080d2917SAlp Dener         }
916080d2917SAlp Dener       }
917080d2917SAlp Dener     } else {
918080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
919080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
920080d2917SAlp Dener     }
921080d2917SAlp Dener     break;
922080d2917SAlp Dener 
923080d2917SAlp Dener   default:
924080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
925b1c2d0e3SAlp Dener       if (prered < 0.0) {
926080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
927080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
928080d2917SAlp Dener         /*  be rejected! */
929080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
930080d2917SAlp Dener       } else {
931b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
932080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
933080d2917SAlp Dener         } else {
934080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
935080d2917SAlp Dener             kappa = 1.0;
936080d2917SAlp Dener           } else {
937080d2917SAlp Dener             kappa = actred / prered;
938080d2917SAlp Dener           }
939080d2917SAlp Dener 
9409566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
941080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
942080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
943080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
944080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
945080d2917SAlp Dener 
946080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
947080d2917SAlp Dener             /*  Great agreement */
948080d2917SAlp Dener             *accept = PETSC_TRUE;
949080d2917SAlp Dener             if (tau_max < 1.0) {
950080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
951080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
952080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
953080d2917SAlp Dener             } else {
954080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
955080d2917SAlp Dener             }
956080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
957080d2917SAlp Dener             /*  Good agreement */
958080d2917SAlp Dener             *accept = PETSC_TRUE;
959080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
960080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
961080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
962080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
963080d2917SAlp Dener             } else if (tau_max < 1.0) {
964080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
965080d2917SAlp Dener             } else {
966080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
967080d2917SAlp Dener             }
968080d2917SAlp Dener           } else {
969080d2917SAlp Dener             /*  Not good agreement */
970080d2917SAlp Dener             if (tau_min > 1.0) {
971080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
972080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
973080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
974080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
975080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
976080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
977080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
978080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
979080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
980080d2917SAlp Dener             } else {
981080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
982080d2917SAlp Dener             }
983080d2917SAlp Dener           }
984080d2917SAlp Dener         }
985080d2917SAlp Dener       }
986080d2917SAlp Dener     } else {
987080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
988080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
989080d2917SAlp Dener     }
99028017e9fSAlp Dener     break;
991080d2917SAlp Dener   }
992c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
993c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
994fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
995080d2917SAlp Dener   PetscFunctionReturn(0);
996080d2917SAlp Dener }
997080d2917SAlp Dener 
998eb910715SAlp Dener /* ---------------------------------------------------------- */
999df278d8fSAlp Dener 
100062675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
100162675beeSAlp Dener {
100262675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
100362675beeSAlp Dener 
100462675beeSAlp Dener   PetscFunctionBegin;
100562675beeSAlp Dener   switch (stepType) {
100662675beeSAlp Dener   case BNK_NEWTON:
100762675beeSAlp Dener     ++bnk->newt;
100862675beeSAlp Dener     break;
100962675beeSAlp Dener   case BNK_BFGS:
101062675beeSAlp Dener     ++bnk->bfgs;
101162675beeSAlp Dener     break;
101262675beeSAlp Dener   case BNK_SCALED_GRADIENT:
101362675beeSAlp Dener     ++bnk->sgrad;
101462675beeSAlp Dener     break;
101562675beeSAlp Dener   case BNK_GRADIENT:
101662675beeSAlp Dener     ++bnk->grad;
101762675beeSAlp Dener     break;
101862675beeSAlp Dener   default:
101962675beeSAlp Dener     break;
102062675beeSAlp Dener   }
102162675beeSAlp Dener   PetscFunctionReturn(0);
102262675beeSAlp Dener }
102362675beeSAlp Dener 
102462675beeSAlp Dener /* ---------------------------------------------------------- */
102562675beeSAlp Dener 
10269b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1027eb910715SAlp Dener {
1028eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1029e031d6f5SAlp Dener   PetscInt       i;
1030eb910715SAlp Dener 
1031eb910715SAlp Dener   PetscFunctionBegin;
1032c4b75bccSAlp Dener   if (!tao->gradient) {
10339566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&tao->gradient));
1034c4b75bccSAlp Dener   }
1035c4b75bccSAlp Dener   if (!tao->stepdirection) {
10369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&tao->stepdirection));
1037c4b75bccSAlp Dener   }
1038c4b75bccSAlp Dener   if (!bnk->W) {
10399566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->W));
1040c4b75bccSAlp Dener   }
1041c4b75bccSAlp Dener   if (!bnk->Xold) {
10429566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->Xold));
1043c4b75bccSAlp Dener   }
1044c4b75bccSAlp Dener   if (!bnk->Gold) {
10459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->Gold));
1046c4b75bccSAlp Dener   }
1047c4b75bccSAlp Dener   if (!bnk->Xwork) {
10489566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->Xwork));
1049c4b75bccSAlp Dener   }
1050c4b75bccSAlp Dener   if (!bnk->Gwork) {
10519566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->Gwork));
1052c4b75bccSAlp Dener   }
1053c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
10549566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->unprojected_gradient));
1055c4b75bccSAlp Dener   }
1056c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
10579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->unprojected_gradient_old));
1058c4b75bccSAlp Dener   }
1059c4b75bccSAlp Dener   if (!bnk->Diag_min) {
10609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->Diag_min));
1061c4b75bccSAlp Dener   }
1062c4b75bccSAlp Dener   if (!bnk->Diag_max) {
10639566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution,&bnk->Diag_max));
1064c4b75bccSAlp Dener   }
1065e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1066c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1067c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
10689566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old)));
10699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old));
107089da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
10719566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient)));
10729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient));
1073c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
10749566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->Gold)));
10759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->G_old));
1076c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
10779566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(tao->gradient)));
10789566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg->gradient));
1079c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
10809566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(tao->stepdirection)));
10819566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg->stepdirection));
1082c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
10839566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(bnk->bncg, tao->solution));
1084c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
10859566063dSJacob Faibussowitsch     PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its));
10869566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol));
10879566063dSJacob Faibussowitsch     PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin));
10889566063dSJacob Faibussowitsch     PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP));
10899566063dSJacob Faibussowitsch     PetscCall(TaoSetObjective(bnk->bncg, tao->ops->computeobjective, tao->user_objP));
10909566063dSJacob Faibussowitsch     PetscCall(TaoSetGradient(bnk->bncg, NULL, tao->ops->computegradient, tao->user_gradP));
10919566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(bnk->bncg, NULL, tao->ops->computeobjectiveandgradient, tao->user_objgradP));
10929566063dSJacob Faibussowitsch     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg)));
1093c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
10949566063dSJacob Faibussowitsch       PetscCall(TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
10959566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(tao->monitorcontext[i])));
1096e031d6f5SAlp Dener     }
1097e031d6f5SAlp Dener   }
109883c8fe1dSLisandro Dalcin   bnk->X_inactive = NULL;
109983c8fe1dSLisandro Dalcin   bnk->G_inactive = NULL;
110083c8fe1dSLisandro Dalcin   bnk->inactive_work = NULL;
110183c8fe1dSLisandro Dalcin   bnk->active_work = NULL;
110283c8fe1dSLisandro Dalcin   bnk->inactive_idx = NULL;
110383c8fe1dSLisandro Dalcin   bnk->active_idx = NULL;
110483c8fe1dSLisandro Dalcin   bnk->active_lower = NULL;
110583c8fe1dSLisandro Dalcin   bnk->active_upper = NULL;
110683c8fe1dSLisandro Dalcin   bnk->active_fixed = NULL;
110783c8fe1dSLisandro Dalcin   bnk->M = NULL;
110883c8fe1dSLisandro Dalcin   bnk->H_inactive = NULL;
110983c8fe1dSLisandro Dalcin   bnk->Hpre_inactive = NULL;
1110eb910715SAlp Dener   PetscFunctionReturn(0);
1111eb910715SAlp Dener }
1112eb910715SAlp Dener 
1113eb910715SAlp Dener /*------------------------------------------------------------*/
1114df278d8fSAlp Dener 
1115e0ed867bSAlp Dener PetscErrorCode TaoDestroy_BNK(Tao tao)
1116eb910715SAlp Dener {
1117eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
1118eb910715SAlp Dener 
1119eb910715SAlp Dener   PetscFunctionBegin;
11209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->W));
11219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Xold));
11229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Gold));
11239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Xwork));
11249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Gwork));
11259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->unprojected_gradient));
11269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->unprojected_gradient_old));
11279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Diag_min));
11289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Diag_max));
11299566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_lower));
11309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_upper));
11319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_fixed));
11329566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_idx));
11339566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->inactive_idx));
11349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
11359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
11369566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&bnk->bncg));
11379566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
1138eb910715SAlp Dener   PetscFunctionReturn(0);
1139eb910715SAlp Dener }
1140eb910715SAlp Dener 
1141eb910715SAlp Dener /*------------------------------------------------------------*/
1142df278d8fSAlp Dener 
1143e0ed867bSAlp Dener PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1144eb910715SAlp Dener {
1145eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1146eb910715SAlp Dener 
1147eb910715SAlp Dener   PetscFunctionBegin;
1148d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization");
11499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL));
11509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL));
11519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL));
11529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL));
11539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL));
11549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL));
11559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL));
11569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL));
11579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL));
11589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL));
11599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL));
11609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL));
11619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL));
11629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL));
11639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL));
11649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL));
11659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL));
11669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL));
11679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2,NULL));
11689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3,NULL));
11699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4,NULL));
11709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5,NULL));
11719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1,NULL));
11729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2,NULL));
11739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3,NULL));
11749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4,NULL));
11759566063dSJacob Faibussowitsch   PetscCall(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));
11769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)", "", bnk->omega2, &bnk->omega2,NULL));
11779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3,NULL));
11789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)", "", bnk->omega4, &bnk->omega4,NULL));
11799566063dSJacob Faibussowitsch   PetscCall(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));
11809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu1_i", "(developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)", "", bnk->mu1_i, &bnk->mu1_i,NULL));
11819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i,NULL));
11829566063dSJacob Faibussowitsch   PetscCall(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));
11839566063dSJacob Faibussowitsch   PetscCall(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));
11849566063dSJacob Faibussowitsch   PetscCall(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));
11859566063dSJacob Faibussowitsch   PetscCall(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));
11869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i,NULL));
11879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL));
11889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL));
11899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1,NULL));
11909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2,NULL));
11919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3,NULL));
11929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4,NULL));
11939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL));
11949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL));
11959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL));
11969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL));
11979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL));
11989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL));
11999566063dSJacob Faibussowitsch   PetscCall(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));
1200d0609cedSBarry Smith   PetscOptionsHeadEnd();
12018ebe3e4eSStefano Zampini 
12029566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(bnk->bncg,((PetscObject)(tao))->prefix));
12039566063dSJacob Faibussowitsch   PetscCall(TaoAppendOptionsPrefix(bnk->bncg,"tao_bnk_cg_"));
12049566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(bnk->bncg));
12058ebe3e4eSStefano Zampini 
12069566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp,((PetscObject)(tao))->prefix));
12079566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp,"tao_bnk_"));
12089566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
1209eb910715SAlp Dener   PetscFunctionReturn(0);
1210eb910715SAlp Dener }
1211eb910715SAlp Dener 
1212eb910715SAlp Dener /*------------------------------------------------------------*/
1213df278d8fSAlp Dener 
1214e0ed867bSAlp Dener PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1215eb910715SAlp Dener {
1216eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1217eb910715SAlp Dener   PetscInt       nrejects;
1218eb910715SAlp Dener   PetscBool      isascii;
1219eb910715SAlp Dener 
1220eb910715SAlp Dener   PetscFunctionBegin;
12219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1222eb910715SAlp Dener   if (isascii) {
12239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1224b9ac7092SAlp Dener     if (bnk->M) {
12259566063dSJacob Faibussowitsch       PetscCall(MatLMVMGetRejectCount(bnk->M,&nrejects));
122663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n",nrejects));
1227eb910715SAlp Dener     }
122863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its));
122963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt));
1230b9ac7092SAlp Dener     if (bnk->M) {
123163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs));
1232b9ac7092SAlp Dener     }
123363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad));
123463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad));
12359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n"));
123663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  atol: %" PetscInt_FMT "\n", bnk->ksp_atol));
123763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol));
123863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol));
123963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  negc: %" PetscInt_FMT "\n", bnk->ksp_negc));
124063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol));
124163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  iter: %" PetscInt_FMT "\n", bnk->ksp_iter));
124263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  othr: %" PetscInt_FMT "\n", bnk->ksp_othr));
12439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1244eb910715SAlp Dener   }
1245eb910715SAlp Dener   PetscFunctionReturn(0);
1246eb910715SAlp Dener }
1247eb910715SAlp Dener 
1248eb910715SAlp Dener /* ---------------------------------------------------------- */
1249df278d8fSAlp Dener 
1250eb910715SAlp Dener /*MC
1251eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
125266ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1253eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1254eb910715SAlp Dener               Hk dk = -gk
12552b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12562b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1257eb910715SAlp Dener 
1258eb910715SAlp Dener     Options Database Keys:
12599fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
12609fa2b5dcSStefano Zampini . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
12619fa2b5dcSStefano Zampini . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
12629fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
12639fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
12649fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
12659fa2b5dcSStefano Zampini . -tao_bnk_sval - (developer) Hessian perturbation starting value
12669fa2b5dcSStefano Zampini . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
12679fa2b5dcSStefano Zampini . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
12689fa2b5dcSStefano Zampini . -tao_bnk_pmin - (developer) minimum Hessian perturbation
12699fa2b5dcSStefano Zampini . -tao_bnk_pmax - (developer) aximum Hessian perturbation
12709fa2b5dcSStefano Zampini . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
12719fa2b5dcSStefano Zampini . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
12729fa2b5dcSStefano Zampini . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
12739fa2b5dcSStefano Zampini . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
12749fa2b5dcSStefano Zampini . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
12759fa2b5dcSStefano Zampini . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction)
12769fa2b5dcSStefano Zampini . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
12779fa2b5dcSStefano Zampini . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
12789fa2b5dcSStefano Zampini . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction)
12799fa2b5dcSStefano Zampini . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
12809fa2b5dcSStefano Zampini . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
12819fa2b5dcSStefano Zampini . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
12829fa2b5dcSStefano Zampini . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
12839fa2b5dcSStefano Zampini . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
12849fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
12859fa2b5dcSStefano Zampini . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
12869fa2b5dcSStefano Zampini . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation)
12879fa2b5dcSStefano Zampini . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
12889fa2b5dcSStefano Zampini . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
12899fa2b5dcSStefano Zampini . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
12909fa2b5dcSStefano Zampini . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
12919fa2b5dcSStefano Zampini . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation)
12929fa2b5dcSStefano Zampini . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step)
12939fa2b5dcSStefano Zampini . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
12949fa2b5dcSStefano Zampini . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step)
12959fa2b5dcSStefano Zampini . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step)
12969fa2b5dcSStefano Zampini . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
12979fa2b5dcSStefano Zampini . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
12989fa2b5dcSStefano Zampini . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step)
12999fa2b5dcSStefano Zampini . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
13009fa2b5dcSStefano Zampini . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
13019fa2b5dcSStefano Zampini . -tao_bnk_mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
13029fa2b5dcSStefano Zampini . -tao_bnk_mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
13039fa2b5dcSStefano Zampini . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
13049fa2b5dcSStefano Zampini . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
13059fa2b5dcSStefano Zampini . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
13069fa2b5dcSStefano Zampini . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
13079fa2b5dcSStefano Zampini - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1308eb910715SAlp Dener 
1309eb910715SAlp Dener   Level: beginner
1310eb910715SAlp Dener M*/
1311eb910715SAlp Dener 
1312eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1313eb910715SAlp Dener {
1314eb910715SAlp Dener   TAO_BNK *bnk;
1315b9ac7092SAlp Dener   PC      pc;
1316eb910715SAlp Dener 
1317eb910715SAlp Dener   PetscFunctionBegin;
13189566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao,&bnk));
1319eb910715SAlp Dener 
1320eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1321eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1322eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1323eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1324eb910715SAlp Dener 
1325eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1326eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1327eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1328eb910715SAlp Dener 
1329eb910715SAlp Dener   tao->data = (void*)bnk;
1330eb910715SAlp Dener 
133166ed3702SAlp Dener   /*  Hessian shifting parameters */
1332e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1333e0ed867bSAlp Dener   bnk->computestep = TaoBNKComputeStep;
1334e0ed867bSAlp Dener 
1335eb910715SAlp Dener   bnk->sval   = 0.0;
1336eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1337eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1338eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1339eb910715SAlp Dener 
1340eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1341eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1342eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1343eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1344eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1345eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1346eb910715SAlp Dener 
1347eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1348eb910715SAlp Dener   bnk->nu1 = 0.25;
1349eb910715SAlp Dener   bnk->nu2 = 0.50;
1350eb910715SAlp Dener   bnk->nu3 = 1.00;
1351eb910715SAlp Dener   bnk->nu4 = 1.25;
1352eb910715SAlp Dener 
1353eb910715SAlp Dener   bnk->omega1 = 0.25;
1354eb910715SAlp Dener   bnk->omega2 = 0.50;
1355eb910715SAlp Dener   bnk->omega3 = 1.00;
1356eb910715SAlp Dener   bnk->omega4 = 2.00;
1357eb910715SAlp Dener   bnk->omega5 = 4.00;
1358eb910715SAlp Dener 
1359eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1360eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1361eb910715SAlp Dener   bnk->eta2 = 0.25;
1362eb910715SAlp Dener   bnk->eta3 = 0.50;
1363eb910715SAlp Dener   bnk->eta4 = 0.90;
1364eb910715SAlp Dener 
1365eb910715SAlp Dener   bnk->alpha1 = 0.25;
1366eb910715SAlp Dener   bnk->alpha2 = 0.50;
1367eb910715SAlp Dener   bnk->alpha3 = 1.00;
1368eb910715SAlp Dener   bnk->alpha4 = 2.00;
1369eb910715SAlp Dener   bnk->alpha5 = 4.00;
1370eb910715SAlp Dener 
1371eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1372eb910715SAlp Dener   bnk->mu1 = 0.10;
1373eb910715SAlp Dener   bnk->mu2 = 0.50;
1374eb910715SAlp Dener 
1375eb910715SAlp Dener   bnk->gamma1 = 0.25;
1376eb910715SAlp Dener   bnk->gamma2 = 0.50;
1377eb910715SAlp Dener   bnk->gamma3 = 2.00;
1378eb910715SAlp Dener   bnk->gamma4 = 4.00;
1379eb910715SAlp Dener 
1380eb910715SAlp Dener   bnk->theta = 0.05;
1381eb910715SAlp Dener 
1382eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1383eb910715SAlp Dener   bnk->mu1_i = 0.35;
1384eb910715SAlp Dener   bnk->mu2_i = 0.50;
1385eb910715SAlp Dener 
1386eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1387eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1388eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1389eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1390eb910715SAlp Dener 
1391eb910715SAlp Dener   bnk->theta_i = 0.25;
1392eb910715SAlp Dener 
1393eb910715SAlp Dener   /*  Remaining parameters */
1394c0f10754SAlp Dener   bnk->max_cg_its = 0;
1395eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1396eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1397770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
13980a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
13990a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
140062675beeSAlp Dener   bnk->dmin = 1.0e-6;
140162675beeSAlp Dener   bnk->dmax = 1.0e6;
1402eb910715SAlp Dener 
140383c8fe1dSLisandro Dalcin   bnk->M           = NULL;
140483c8fe1dSLisandro Dalcin   bnk->bfgs_pre    = NULL;
1405eb910715SAlp Dener   bnk->init_type   = BNK_INIT_INTERPOLATION;
14067b1c7716SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION;
14072f75a4aaSAlp Dener   bnk->as_type     = BNK_AS_BERTSEKAS;
1408eb910715SAlp Dener 
1409e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
14109566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&bnk->bncg));
14119566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg,(PetscObject)tao,1));
14129566063dSJacob Faibussowitsch   PetscCall(TaoSetType(bnk->bncg,TAOBNCG));
1413e031d6f5SAlp Dener 
1414c0f10754SAlp Dener   /* Create the line search */
14159566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch));
14169566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch,(PetscObject)tao,1));
1417*f4db9bf7SStefano Zampini   PetscCall(TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT));
14189566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao));
1419eb910715SAlp Dener 
1420eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
14219566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp));
14229566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp,(PetscObject)tao,1));
14239566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp,KSPSTCG));
14249566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp,&pc));
14259566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCLMVM));
1426eb910715SAlp Dener   PetscFunctionReturn(0);
1427eb910715SAlp Dener }
1428