xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 976ed0a4a7e2c9a504a034b07aaf489e2c2d55c5)
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 
9b3e6a353SBarry Smith /* Extracts from the full Hessian the part associated with the current bnk->inactive_idx and set the PCLMVM preconditioner */
10e031d6f5SAlp Dener 
11b3e6a353SBarry Smith static PetscErrorCode TaoBNKComputeSubHessian(Tao tao)
12b3e6a353SBarry Smith {
13b3e6a353SBarry Smith   TAO_BNK *bnk = (TAO_BNK *)tao->data;
14b3e6a353SBarry Smith 
15b3e6a353SBarry Smith   PetscFunctionBegin;
16b3e6a353SBarry Smith   PetscCall(MatDestroy(&bnk->Hpre_inactive));
17b3e6a353SBarry Smith   PetscCall(MatDestroy(&bnk->H_inactive));
18b3e6a353SBarry Smith   if (bnk->active_idx) {
19b3e6a353SBarry Smith     PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
20b3e6a353SBarry Smith     if (tao->hessian == tao->hessian_pre) {
21b3e6a353SBarry Smith       PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
22b3e6a353SBarry Smith       bnk->Hpre_inactive = bnk->H_inactive;
23b3e6a353SBarry Smith     } else {
24b3e6a353SBarry Smith       PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive));
25b3e6a353SBarry Smith     }
26b3e6a353SBarry Smith     if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx));
27b3e6a353SBarry Smith   } else {
28b3e6a353SBarry Smith     PetscCall(PetscObjectReference((PetscObject)tao->hessian));
29b3e6a353SBarry Smith     bnk->H_inactive = tao->hessian;
30b3e6a353SBarry Smith     PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre));
31b3e6a353SBarry Smith     bnk->Hpre_inactive = tao->hessian_pre;
32b3e6a353SBarry Smith     if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre));
33b3e6a353SBarry Smith   }
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35b3e6a353SBarry Smith }
36b3e6a353SBarry Smith 
37b3e6a353SBarry Smith /* Initializes the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
38df278d8fSAlp Dener 
39d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
40d71ae5a4SJacob Faibussowitsch {
41eb910715SAlp Dener   TAO_BNK          *bnk = (TAO_BNK *)tao->data;
42eb910715SAlp Dener   PC                pc;
4389da521bSAlp Dener   PetscReal         f_min, ftrial, prered, actred, kappa, sigma, resnorm;
44eb910715SAlp Dener   PetscReal         tau, tau_1, tau_2, tau_max, tau_min, max_radius;
450ad3a497SAlp Dener   PetscBool         is_bfgs, is_jacobi, is_symmetric, sym_set;
46c4b75bccSAlp Dener   PetscInt          n, N, nDiff;
47eb910715SAlp Dener   PetscInt          i_max = 5;
48eb910715SAlp Dener   PetscInt          j_max = 1;
49eb910715SAlp Dener   PetscInt          i, j;
502e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
51eb910715SAlp Dener 
52eb910715SAlp Dener   PetscFunctionBegin;
5328017e9fSAlp Dener   /* Project the current point onto the feasible set */
549566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
559566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU));
561baa6e33SBarry Smith   if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
5728017e9fSAlp Dener 
5828017e9fSAlp Dener   /* Project the initial point onto the feasible region */
599566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
6028017e9fSAlp Dener 
6128017e9fSAlp Dener   /* Check convergence criteria */
629566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
639566063dSJacob Faibussowitsch   PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
649566063dSJacob Faibussowitsch   PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
65*976ed0a4SStefano Zampini   if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
669566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
6728017e9fSAlp Dener 
68c0f10754SAlp Dener   /* Test the initial point for convergence */
699566063dSJacob Faibussowitsch   PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
709566063dSJacob Faibussowitsch   PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
713c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
729566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
739566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
74dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
753ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
76c0f10754SAlp Dener 
77e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
78eb910715SAlp Dener   bnk->ksp_atol = 0;
79eb910715SAlp Dener   bnk->ksp_rtol = 0;
80eb910715SAlp Dener   bnk->ksp_dtol = 0;
81eb910715SAlp Dener   bnk->ksp_ctol = 0;
82eb910715SAlp Dener   bnk->ksp_negc = 0;
83eb910715SAlp Dener   bnk->ksp_iter = 0;
84eb910715SAlp Dener   bnk->ksp_othr = 0;
85eb910715SAlp Dener 
86e031d6f5SAlp Dener   /* Reset accepted step type counters */
87e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
88e031d6f5SAlp Dener   bnk->newt       = 0;
89e031d6f5SAlp Dener   bnk->bfgs       = 0;
90e031d6f5SAlp Dener   bnk->sgrad      = 0;
91e031d6f5SAlp Dener   bnk->grad       = 0;
92e031d6f5SAlp Dener 
93fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
94fed79b8eSAlp Dener   bnk->pert = bnk->sval;
95fed79b8eSAlp Dener 
96937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
979566063dSJacob Faibussowitsch   PetscCall(VecSet(tao->stepdirection, 0.0));
98937a31a1SAlp Dener 
99e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
1009566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
103b9ac7092SAlp Dener   if (is_bfgs) {
104b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
1059566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M));
1069566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
1079566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
1089566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(bnk->M, n, n, N, N));
1099566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient));
1109566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric));
1113c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
1121baa6e33SBarry Smith   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
113e031d6f5SAlp Dener 
114e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
1159566063dSJacob Faibussowitsch   PetscCall(VecSet(bnk->Diag_min, bnk->dmin));
1169566063dSJacob Faibussowitsch   PetscCall(VecSet(bnk->Diag_max, bnk->dmax));
117eb910715SAlp Dener 
118eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
119eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
120c0f10754SAlp Dener   *needH = PETSC_TRUE;
1219566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp, "KSPCGSetRadius_C", &kspTR));
1222e6e4ca1SStefano Zampini   if (kspTR) {
12362675beeSAlp Dener     switch (initType) {
124eb910715SAlp Dener     case BNK_INIT_CONSTANT:
125eb910715SAlp Dener       /* Use the initial radius specified */
126c0f10754SAlp Dener       tao->trust = tao->trust0;
127eb910715SAlp Dener       break;
128eb910715SAlp Dener 
129eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
130c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
131eb910715SAlp Dener       max_radius = 0.0;
13208752603SAlp Dener       tao->trust = tao->trust0;
133eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1340a4511e9SAlp Dener         f_min = bnk->f;
135eb910715SAlp Dener         sigma = 0.0;
136eb910715SAlp Dener 
137c0f10754SAlp Dener         if (*needH) {
13862602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
1399566063dSJacob Faibussowitsch           PetscCall((*bnk->computehessian)(tao));
1409566063dSJacob Faibussowitsch           PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE));
141b3e6a353SBarry Smith           PetscCall(TaoBNKComputeSubHessian(tao));
142c0f10754SAlp Dener           *needH = PETSC_FALSE;
143eb910715SAlp Dener         }
144eb910715SAlp Dener 
145eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
14662602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
1479566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->Xold));
1489566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, -tao->trust / bnk->gnorm, tao->gradient));
1499566063dSJacob Faibussowitsch           PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
15089da521bSAlp Dener           /* Compute the step we actually accepted */
1519566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->W));
1529566063dSJacob Faibussowitsch           PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold));
15362602cfbSAlp Dener           /* Compute the objective at the trial */
1549566063dSJacob Faibussowitsch           PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial));
1553c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
1569566063dSJacob Faibussowitsch           PetscCall(VecCopy(bnk->Xold, tao->solution));
157eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
158eb910715SAlp Dener             tau = bnk->gamma1_i;
159eb910715SAlp Dener           } else {
1600a4511e9SAlp Dener             if (ftrial < f_min) {
1610a4511e9SAlp Dener               f_min = ftrial;
162eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
163eb910715SAlp Dener             }
16408752603SAlp Dener 
165770b7498SAlp Dener             /* Compute the predicted and actual reduction */
16689da521bSAlp Dener             if (bnk->active_idx) {
1679566063dSJacob Faibussowitsch               PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1689566063dSJacob Faibussowitsch               PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1692ab2a32cSAlp Dener             } else {
17008752603SAlp Dener               bnk->X_inactive    = bnk->W;
17108752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1722ab2a32cSAlp Dener             }
1739566063dSJacob Faibussowitsch             PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
1749566063dSJacob Faibussowitsch             PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered));
17589da521bSAlp Dener             if (bnk->active_idx) {
1769566063dSJacob Faibussowitsch               PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1779566063dSJacob Faibussowitsch               PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1782ab2a32cSAlp Dener             }
179eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
180eb910715SAlp Dener             actred = bnk->f - ftrial;
1813105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
182eb910715SAlp Dener               kappa = 1.0;
1833105154fSTodd Munson             } else {
184eb910715SAlp Dener               kappa = actred / prered;
185eb910715SAlp Dener             }
186eb910715SAlp Dener 
187eb910715SAlp Dener             tau_1   = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
188eb910715SAlp Dener             tau_2   = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
189eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
190eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
191eb910715SAlp Dener 
19218cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
193eb910715SAlp Dener               /*  Great agreement */
194eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
195eb910715SAlp Dener 
196eb910715SAlp Dener               if (tau_max < 1.0) {
197eb910715SAlp Dener                 tau = bnk->gamma3_i;
1983105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
199eb910715SAlp Dener                 tau = bnk->gamma4_i;
2003105154fSTodd Munson               } else {
201eb910715SAlp Dener                 tau = tau_max;
202eb910715SAlp Dener               }
20318cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
204eb910715SAlp Dener               /*  Good agreement */
205eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
206eb910715SAlp Dener 
207eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
208eb910715SAlp Dener                 tau = bnk->gamma2_i;
209eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
210eb910715SAlp Dener                 tau = bnk->gamma3_i;
211eb910715SAlp Dener               } else {
212eb910715SAlp Dener                 tau = tau_max;
213eb910715SAlp Dener               }
2148f8a4e06SAlp Dener             } else {
215eb910715SAlp Dener               /*  Not good agreement */
216eb910715SAlp Dener               if (tau_min > 1.0) {
217eb910715SAlp Dener                 tau = bnk->gamma2_i;
218eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
219eb910715SAlp Dener                 tau = bnk->gamma1_i;
220eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
221eb910715SAlp Dener                 tau = bnk->gamma1_i;
2223105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
223eb910715SAlp Dener                 tau = tau_1;
2243105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
225eb910715SAlp Dener                 tau = tau_2;
226eb910715SAlp Dener               } else {
227eb910715SAlp Dener                 tau = tau_max;
228eb910715SAlp Dener               }
229eb910715SAlp Dener             }
230eb910715SAlp Dener           }
231eb910715SAlp Dener           tao->trust = tau * tao->trust;
232eb910715SAlp Dener         }
233eb910715SAlp Dener 
2340a4511e9SAlp Dener         if (f_min < bnk->f) {
235937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2360a4511e9SAlp Dener           bnk->f = f_min;
2379566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->Xold));
2389566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
2399566063dSJacob Faibussowitsch           PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
2409566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, tao->stepdirection));
2419566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
2429566063dSJacob Faibussowitsch           PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
2439566063dSJacob Faibussowitsch           PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
2449566063dSJacob Faibussowitsch           PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
245*976ed0a4SStefano Zampini           if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
246937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
2479566063dSJacob Faibussowitsch           PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
248c0f10754SAlp Dener           *needH = PETSC_TRUE;
249937a31a1SAlp Dener           /* Test the new step for convergence */
2509566063dSJacob Faibussowitsch           PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
2519566063dSJacob Faibussowitsch           PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
2523c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2539566063dSJacob Faibussowitsch           PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
2549566063dSJacob Faibussowitsch           PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
255dbbe0bcdSBarry Smith           PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
2563ba16761SJacob Faibussowitsch           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
257937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
2589566063dSJacob Faibussowitsch           PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
259eb910715SAlp Dener         }
260eb910715SAlp Dener       }
261eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
262e031d6f5SAlp Dener 
263e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
264e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
265e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
266eb910715SAlp Dener       break;
267eb910715SAlp Dener 
268eb910715SAlp Dener     default:
269eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
270eb910715SAlp Dener       tao->trust = 0.0;
271eb910715SAlp Dener       break;
272eb910715SAlp Dener     }
273eb910715SAlp Dener   }
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275eb910715SAlp Dener }
276eb910715SAlp Dener 
277df278d8fSAlp Dener /*------------------------------------------------------------*/
278df278d8fSAlp Dener 
279b3e6a353SBarry Smith /* Computes the exact Hessian and extracts its subHessian */
28062675beeSAlp Dener 
281d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKComputeHessian(Tao tao)
282d71ae5a4SJacob Faibussowitsch {
28362675beeSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
28462675beeSAlp Dener 
28562675beeSAlp Dener   PetscFunctionBegin;
28662675beeSAlp Dener   /* Compute the Hessian */
2879566063dSJacob Faibussowitsch   PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
28862675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
2891baa6e33SBarry Smith   if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
290e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
291b3e6a353SBarry Smith   PetscCall(TaoBNKComputeSubHessian(tao));
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29362675beeSAlp Dener }
29462675beeSAlp Dener 
29562675beeSAlp Dener /*------------------------------------------------------------*/
29662675beeSAlp Dener 
2972f75a4aaSAlp Dener /* Routine for estimating the active set */
2982f75a4aaSAlp Dener 
299d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
300d71ae5a4SJacob Faibussowitsch {
3012f75a4aaSAlp Dener   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
302f4db9bf7SStefano Zampini   PetscBool hessComputed, diagExists, hadactive;
3032f75a4aaSAlp Dener 
3042f75a4aaSAlp Dener   PetscFunctionBegin;
305f4db9bf7SStefano Zampini   hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE;
30608752603SAlp Dener   switch (asType) {
3072f75a4aaSAlp Dener   case BNK_AS_NONE:
3089566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bnk->inactive_idx));
3099566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx));
3109566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bnk->active_idx));
3119566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx));
3122f75a4aaSAlp Dener     break;
3132f75a4aaSAlp Dener 
3142f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3152f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
316b9ac7092SAlp Dener     if (bnk->M) {
3172f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3189566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W));
3192f75a4aaSAlp Dener     } else {
320fc5ca067SStefano Zampini       hessComputed = diagExists = PETSC_FALSE;
32148a46eb9SPierre Jolivet       if (tao->hessian) PetscCall(MatAssembled(tao->hessian, &hessComputed));
32248a46eb9SPierre Jolivet       if (hessComputed) PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists));
323fc5ca067SStefano Zampini       if (diagExists) {
3249b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3259566063dSJacob Faibussowitsch         PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork));
3269566063dSJacob Faibussowitsch         PetscCall(VecAbs(bnk->Xwork));
3279566063dSJacob Faibussowitsch         PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork));
3289566063dSJacob Faibussowitsch         PetscCall(VecReciprocal(bnk->Xwork));
3299566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient));
33061be54a6SAlp Dener       } else {
331c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
3329566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W));
33361be54a6SAlp Dener       }
3342f75a4aaSAlp Dener     }
3359566063dSJacob Faibussowitsch     PetscCall(VecScale(bnk->W, -1.0));
3369371c9d4SSatish Balay     PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx));
337c4b75bccSAlp Dener     break;
3382f75a4aaSAlp Dener 
339d71ae5a4SJacob Faibussowitsch   default:
340d71ae5a4SJacob Faibussowitsch     break;
3412f75a4aaSAlp Dener   }
342f4db9bf7SStefano Zampini   bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3442f75a4aaSAlp Dener }
3452f75a4aaSAlp Dener 
3462f75a4aaSAlp Dener /*------------------------------------------------------------*/
3472f75a4aaSAlp Dener 
3482f75a4aaSAlp Dener /* Routine for bounding the step direction */
3492f75a4aaSAlp Dener 
350d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
351d71ae5a4SJacob Faibussowitsch {
3522f75a4aaSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
3532f75a4aaSAlp Dener 
3542f75a4aaSAlp Dener   PetscFunctionBegin;
355a1318120SAlp Dener   switch (asType) {
356d71ae5a4SJacob Faibussowitsch   case BNK_AS_NONE:
357*976ed0a4SStefano Zampini     if (bnk->active_idx) PetscCall(VecISSet(step, bnk->active_idx, 0.0));
358d71ae5a4SJacob Faibussowitsch     break;
359d71ae5a4SJacob Faibussowitsch   case BNK_AS_BERTSEKAS:
360d71ae5a4SJacob Faibussowitsch     PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step));
361d71ae5a4SJacob Faibussowitsch     break;
362d71ae5a4SJacob Faibussowitsch   default:
363d71ae5a4SJacob Faibussowitsch     break;
3642f75a4aaSAlp Dener   }
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3662f75a4aaSAlp Dener }
3672f75a4aaSAlp Dener 
368e031d6f5SAlp Dener /*------------------------------------------------------------*/
369e031d6f5SAlp Dener 
370e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
371e031d6f5SAlp Dener    accelerate Newton convergence.
372e031d6f5SAlp Dener 
373e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
374e031d6f5SAlp Dener    for more gradient evaluations.
375e031d6f5SAlp Dener */
376e031d6f5SAlp Dener 
377d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
378d71ae5a4SJacob Faibussowitsch {
379c0f10754SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
380c0f10754SAlp Dener 
381c0f10754SAlp Dener   PetscFunctionBegin;
382c0f10754SAlp Dener   *terminate = PETSC_FALSE;
383c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
384c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
385c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
386c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
3879566063dSJacob Faibussowitsch     PetscCall(TaoSolve(bnk->bncg));
388c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
389c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
390c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
391c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
392c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
393e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
394c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
395c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
396c0f10754SAlp 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) {
397c0f10754SAlp Dener       *terminate = PETSC_TRUE;
39861be54a6SAlp Dener     } else {
3999566063dSJacob Faibussowitsch       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
400c0f10754SAlp Dener     }
401c0f10754SAlp Dener   }
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
403c0f10754SAlp Dener }
404c0f10754SAlp Dener 
4052f75a4aaSAlp Dener /*------------------------------------------------------------*/
4062f75a4aaSAlp Dener 
407c0f10754SAlp Dener /* Routine for computing the Newton step. */
408df278d8fSAlp Dener 
409d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
410d71ae5a4SJacob Faibussowitsch {
411eb910715SAlp Dener   TAO_BNK          *bnk         = (TAO_BNK *)tao->data;
412eb910715SAlp Dener   PetscInt          bfgsUpdates = 0;
413eb910715SAlp Dener   PetscInt          kspits;
414bddd1ffdSAlp Dener   PetscBool         is_lmvm;
4152e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
416eb910715SAlp Dener 
417eb910715SAlp Dener   PetscFunctionBegin;
41889da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
41989da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
42089da521bSAlp Dener   if (!bnk->inactive_idx) {
4219566063dSJacob Faibussowitsch     PetscCall(VecSet(tao->stepdirection, 0.0));
4229566063dSJacob Faibussowitsch     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
4233ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
42489da521bSAlp Dener   }
42589da521bSAlp Dener 
42662675beeSAlp Dener   /* Shift the reduced Hessian matrix */
427e831869dSStefano Zampini   if (shift && bnk->pert > 0) {
4289566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm));
429f7bf01afSAlp Dener     if (is_lmvm) {
4309566063dSJacob Faibussowitsch       PetscCall(MatShift(tao->hessian, bnk->pert));
431f7bf01afSAlp Dener     } else {
4329566063dSJacob Faibussowitsch       PetscCall(MatShift(bnk->H_inactive, bnk->pert));
43348a46eb9SPierre Jolivet       if (bnk->H_inactive != bnk->Hpre_inactive) PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert));
43462675beeSAlp Dener     }
435f7bf01afSAlp Dener   }
43662675beeSAlp Dener 
437eb910715SAlp Dener   /* Solve the Newton system of equations */
438937a31a1SAlp Dener   tao->ksp_its = 0;
4399566063dSJacob Faibussowitsch   PetscCall(VecSet(tao->stepdirection, 0.0));
440f4db9bf7SStefano Zampini   if (bnk->resetksp) {
4419566063dSJacob Faibussowitsch     PetscCall(KSPReset(tao->ksp));
4429566063dSJacob Faibussowitsch     PetscCall(KSPResetFromOptions(tao->ksp));
443f4db9bf7SStefano Zampini     bnk->resetksp = PETSC_FALSE;
444f4db9bf7SStefano Zampini   }
4459566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(tao->ksp, bnk->H_inactive, bnk->Hpre_inactive));
4469566063dSJacob Faibussowitsch   PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork));
44789da521bSAlp Dener   if (bnk->active_idx) {
4489566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4499566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4505e9b73cbSAlp Dener   } else {
4515e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4525e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
45328017e9fSAlp Dener   }
4549566063dSJacob Faibussowitsch   PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
4559566063dSJacob Faibussowitsch   PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4569566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
457eb910715SAlp Dener   tao->ksp_its += kspits;
458eb910715SAlp Dener   tao->ksp_tot_its += kspits;
459f4db9bf7SStefano Zampini   PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp, "KSPCGGetNormD_C", &kspTR));
460f4db9bf7SStefano Zampini   if (kspTR) {
4619566063dSJacob Faibussowitsch     PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
462eb910715SAlp Dener 
463eb910715SAlp Dener     if (0.0 == tao->trust) {
464eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
465080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
466080d2917SAlp Dener         tao->trust = bnk->dnorm;
467eb910715SAlp Dener 
468eb910715SAlp Dener         /* Modify the radius if it is too large or small */
469eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
470eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
471eb910715SAlp Dener       } else {
472eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
473eb910715SAlp Dener            the trust-region subproblem to get a direction */
474eb910715SAlp Dener         tao->trust = tao->trust0;
475eb910715SAlp Dener 
476eb910715SAlp Dener         /* Modify the radius if it is too large or small */
477eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
478eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
479eb910715SAlp Dener 
4809566063dSJacob Faibussowitsch         PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
4819566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4829566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
483eb910715SAlp Dener         tao->ksp_its += kspits;
484eb910715SAlp Dener         tao->ksp_tot_its += kspits;
4859566063dSJacob Faibussowitsch         PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
486eb910715SAlp Dener 
4873c859ba3SBarry Smith         PetscCheck(bnk->dnorm != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
488eb910715SAlp Dener       }
489eb910715SAlp Dener     }
490eb910715SAlp Dener   }
4915e9b73cbSAlp Dener   /* Restore sub vectors back */
49289da521bSAlp Dener   if (bnk->active_idx) {
4939566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4949566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4955e9b73cbSAlp Dener   }
496770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
4979566063dSJacob Faibussowitsch   PetscCall(VecScale(tao->stepdirection, -1.0));
4989566063dSJacob Faibussowitsch   PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
499770b7498SAlp Dener 
500770b7498SAlp Dener   /* Record convergence reasons */
5019566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason));
502e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
503770b7498SAlp Dener     ++bnk->ksp_atol;
504e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
505770b7498SAlp Dener     ++bnk->ksp_rtol;
5064a221d59SStefano Zampini   } else if (KSP_CONVERGED_STEP_LENGTH == *ksp_reason) {
507770b7498SAlp Dener     ++bnk->ksp_ctol;
5084a221d59SStefano Zampini   } else if (KSP_CONVERGED_NEG_CURVE == *ksp_reason) {
509770b7498SAlp Dener     ++bnk->ksp_negc;
510e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
511770b7498SAlp Dener     ++bnk->ksp_dtol;
512e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
513770b7498SAlp Dener     ++bnk->ksp_iter;
514770b7498SAlp Dener   } else {
515770b7498SAlp Dener     ++bnk->ksp_othr;
516770b7498SAlp Dener   }
517fed79b8eSAlp Dener 
518fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
519b9ac7092SAlp Dener   if (bnk->M) {
5209566063dSJacob Faibussowitsch     PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
521b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
522fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
5239566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
5249566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
525eb910715SAlp Dener     }
526fed79b8eSAlp Dener   }
5276b591159SAlp Dener   *step_type = BNK_NEWTON;
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
529e465cd6fSAlp Dener }
530eb910715SAlp Dener 
53162675beeSAlp Dener /*------------------------------------------------------------*/
53262675beeSAlp Dener 
5335e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5345e9b73cbSAlp Dener 
535d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
536d71ae5a4SJacob Faibussowitsch {
5375e9b73cbSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
5385e9b73cbSAlp Dener 
5395e9b73cbSAlp Dener   PetscFunctionBegin;
5405e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
54189da521bSAlp Dener   if (bnk->active_idx) {
5429566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5439566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5449566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5455e9b73cbSAlp Dener   } else {
5465e9b73cbSAlp Dener     bnk->X_inactive    = tao->stepdirection;
5475e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5485e9b73cbSAlp Dener     bnk->G_inactive    = bnk->Gwork;
5495e9b73cbSAlp Dener   }
5505e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5519566063dSJacob Faibussowitsch   PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
5529566063dSJacob Faibussowitsch   PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive));
5539566063dSJacob Faibussowitsch   PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered));
5545e9b73cbSAlp Dener   /* Restore the sub vectors */
55589da521bSAlp Dener   if (bnk->active_idx) {
5569566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5579566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5589566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5595e9b73cbSAlp Dener   }
5603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5615e9b73cbSAlp Dener }
5625e9b73cbSAlp Dener 
5635e9b73cbSAlp Dener /*------------------------------------------------------------*/
5645e9b73cbSAlp Dener 
56562675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
56662675beeSAlp Dener 
56762675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
56862675beeSAlp Dener    in the event that the Newton step fails the test.
56962675beeSAlp Dener */
57062675beeSAlp Dener 
571d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
572d71ae5a4SJacob Faibussowitsch {
573e465cd6fSAlp Dener   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
574b2d8c577SAlp Dener   PetscReal gdx, e_min;
575e465cd6fSAlp Dener   PetscInt  bfgsUpdates;
576e465cd6fSAlp Dener 
577e465cd6fSAlp Dener   PetscFunctionBegin;
5786b591159SAlp Dener   switch (*stepType) {
5796b591159SAlp Dener   case BNK_NEWTON:
5809566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
581eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
582eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
583eb910715SAlp Dener         Update the perturbation for next time */
584eb910715SAlp Dener       if (bnk->pert <= 0.0) {
5852e6e4ca1SStefano Zampini         PetscBool is_gltr;
5862e6e4ca1SStefano Zampini 
587eb910715SAlp Dener         /* Initialize the perturbation */
588eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
5899566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr));
5902e6e4ca1SStefano Zampini         if (is_gltr) {
5919566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
592eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
593eb910715SAlp Dener         }
594eb910715SAlp Dener       } else {
595eb910715SAlp Dener         /* Increase the perturbation */
596eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
597eb910715SAlp Dener       }
598eb910715SAlp Dener 
5990ad3a497SAlp Dener       if (!bnk->M) {
600eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
601eb910715SAlp Dener           Must use gradient direction in this case */
6029566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->gradient, tao->stepdirection));
603eb910715SAlp Dener         *stepType = BNK_GRADIENT;
604eb910715SAlp Dener       } else {
605eb910715SAlp Dener         /* Attempt to use the BFGS direction */
6069566063dSJacob Faibussowitsch         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
607eb910715SAlp Dener 
6088d5ead36SAlp Dener         /* Check for success (descent direction)
6098d5ead36SAlp Dener           NOTE: Negative gdx here means not a descent direction because
6108d5ead36SAlp Dener           the fall-back step is missing a negative sign. */
6119566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
6123105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
613eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
614eb910715SAlp Dener             We can assert bfgsUpdates > 1 in this case because
615eb910715SAlp Dener             the first solve produces the scaled gradient direction,
616eb910715SAlp Dener             which is guaranteed to be descent */
617eb910715SAlp Dener 
618eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
6199566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6209566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6219566063dSJacob Faibussowitsch           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
622eb910715SAlp Dener 
623eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
624eb910715SAlp Dener         } else {
6259566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
626eb910715SAlp Dener           if (1 == bfgsUpdates) {
627eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
628eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
629eb910715SAlp Dener           } else {
630eb910715SAlp Dener             *stepType = BNK_BFGS;
631eb910715SAlp Dener           }
632eb910715SAlp Dener         }
633eb910715SAlp Dener       }
6348d5ead36SAlp Dener       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6359566063dSJacob Faibussowitsch       PetscCall(VecScale(tao->stepdirection, -1.0));
6369566063dSJacob Faibussowitsch       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
637eb910715SAlp Dener     } else {
638eb910715SAlp Dener       /* Computed Newton step is descent */
639eb910715SAlp Dener       switch (ksp_reason) {
640eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
641eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
642eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
643eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
6444a221d59SStefano Zampini       case KSP_CONVERGED_NEG_CURVE:
645eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
646eb910715SAlp Dener         if (bnk->pert <= 0.0) {
6472e6e4ca1SStefano Zampini           PetscBool is_gltr;
6482e6e4ca1SStefano Zampini 
649eb910715SAlp Dener           /* Initialize the perturbation */
650eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
6519566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr));
6522e6e4ca1SStefano Zampini           if (is_gltr) {
6539566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
654eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
655eb910715SAlp Dener           }
656eb910715SAlp Dener         } else {
657eb910715SAlp Dener           /* Increase the perturbation */
658eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
659eb910715SAlp Dener         }
660eb910715SAlp Dener         break;
661eb910715SAlp Dener 
662eb910715SAlp Dener       default:
663eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
664eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
665ad540459SPierre Jolivet         if (bnk->pert < bnk->pmin) bnk->pert = 0.0;
666eb910715SAlp Dener         break;
667eb910715SAlp Dener       }
668fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
669eb910715SAlp Dener     }
6706b591159SAlp Dener     break;
6716b591159SAlp Dener 
6726b591159SAlp Dener   case BNK_BFGS:
6736b591159SAlp Dener     /* Check for success (descent direction) */
6749566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
6756b591159SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
6766b591159SAlp Dener       /* Step is not descent or solve was not successful
6776b591159SAlp Dener          Use steepest descent direction (scaled) */
6789566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6799566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6809566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection));
6819566063dSJacob Faibussowitsch       PetscCall(VecScale(tao->stepdirection, -1.0));
6829566063dSJacob Faibussowitsch       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
6836b591159SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
6846b591159SAlp Dener     } else {
6856b591159SAlp Dener       *stepType = BNK_BFGS;
6866b591159SAlp Dener     }
6876b591159SAlp Dener     break;
6886b591159SAlp Dener 
689d71ae5a4SJacob Faibussowitsch   case BNK_SCALED_GRADIENT:
690d71ae5a4SJacob Faibussowitsch     break;
6916b591159SAlp Dener 
692d71ae5a4SJacob Faibussowitsch   default:
693d71ae5a4SJacob Faibussowitsch     break;
6946b591159SAlp Dener   }
6956b591159SAlp Dener 
6963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
697eb910715SAlp Dener }
698eb910715SAlp Dener 
699df278d8fSAlp Dener /*------------------------------------------------------------*/
700df278d8fSAlp Dener 
701df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
702df278d8fSAlp Dener 
703df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
704df278d8fSAlp Dener   Newton step does not produce a valid step length.
705df278d8fSAlp Dener */
706df278d8fSAlp Dener 
707d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
708d71ae5a4SJacob Faibussowitsch {
709c14b763aSAlp Dener   TAO_BNK                     *bnk = (TAO_BNK *)tao->data;
710c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
711b2d8c577SAlp Dener   PetscReal                    e_min, gdx;
712c14b763aSAlp Dener   PetscInt                     bfgsUpdates;
713c14b763aSAlp Dener 
714c14b763aSAlp Dener   PetscFunctionBegin;
715c14b763aSAlp Dener   /* Perform the linesearch */
7169566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7179566063dSJacob Faibussowitsch   PetscCall(TaoAddLineSearchCounts(tao));
718c14b763aSAlp Dener 
719b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
720c14b763aSAlp Dener     /* Linesearch failed, revert solution */
721c14b763aSAlp Dener     bnk->f = bnk->fold;
7229566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->Xold, tao->solution));
7239566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
724c14b763aSAlp Dener 
725937a31a1SAlp Dener     switch (*stepType) {
726c14b763aSAlp Dener     case BNK_NEWTON:
7278d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
728c14b763aSAlp Dener          Update the perturbation for next time */
729c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
7302e6e4ca1SStefano Zampini         PetscBool is_gltr;
7312e6e4ca1SStefano Zampini 
732c14b763aSAlp Dener         /* Initialize the perturbation */
733c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
7349566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr));
7352e6e4ca1SStefano Zampini         if (is_gltr) {
7369566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
737c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
738c14b763aSAlp Dener         }
739c14b763aSAlp Dener       } else {
740c14b763aSAlp Dener         /* Increase the perturbation */
741c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
742c14b763aSAlp Dener       }
743c14b763aSAlp Dener 
7440ad3a497SAlp Dener       if (!bnk->M) {
745c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
746c14b763aSAlp Dener            Must use gradient direction in this case */
7479566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection));
748937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
749c14b763aSAlp Dener       } else {
750c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
7519566063dSJacob Faibussowitsch         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
7528d5ead36SAlp Dener         /* Check for success (descent direction)
7538d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
7549566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
7553105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
756c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
757c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
758c14b763aSAlp Dener              Use steepest descent direction (scaled) */
7599566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7609566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7619566063dSJacob Faibussowitsch           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
762c14b763aSAlp Dener 
763c14b763aSAlp Dener           bfgsUpdates = 1;
764937a31a1SAlp Dener           *stepType   = BNK_SCALED_GRADIENT;
765c14b763aSAlp Dener         } else {
7669566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
767c14b763aSAlp Dener           if (1 == bfgsUpdates) {
768c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
769937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
770c14b763aSAlp Dener           } else {
771937a31a1SAlp Dener             *stepType = BNK_BFGS;
772c14b763aSAlp Dener           }
773c14b763aSAlp Dener         }
774c14b763aSAlp Dener       }
775c14b763aSAlp Dener       break;
776c14b763aSAlp Dener 
777c14b763aSAlp Dener     case BNK_BFGS:
778c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
779c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
780c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
7819566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7829566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7839566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
784c14b763aSAlp Dener 
785c14b763aSAlp Dener       bfgsUpdates = 1;
786937a31a1SAlp Dener       *stepType   = BNK_SCALED_GRADIENT;
787c14b763aSAlp Dener       break;
788c14b763aSAlp Dener     }
7898d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7909566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
7919566063dSJacob Faibussowitsch     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
792c14b763aSAlp Dener 
7938d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
7949566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7959566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
796c14b763aSAlp Dener   }
797c14b763aSAlp Dener   *reason = ls_reason;
7983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
799c14b763aSAlp Dener }
800c14b763aSAlp Dener 
801df278d8fSAlp Dener /*------------------------------------------------------------*/
802df278d8fSAlp Dener 
803df278d8fSAlp Dener /* Routine for updating the trust radius.
804df278d8fSAlp Dener 
805df278d8fSAlp Dener   Function features three different update methods:
806df278d8fSAlp Dener   1) Line-search step length based
807df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
808df278d8fSAlp Dener   3) Interpolation
809df278d8fSAlp Dener */
810df278d8fSAlp Dener 
811d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
812d71ae5a4SJacob Faibussowitsch {
813080d2917SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
814080d2917SAlp Dener 
815b1c2d0e3SAlp Dener   PetscReal step, kappa;
816080d2917SAlp Dener   PetscReal gdx, tau_1, tau_2, tau_min, tau_max;
817080d2917SAlp Dener 
818080d2917SAlp Dener   PetscFunctionBegin;
819080d2917SAlp Dener   /* Update trust region radius */
820080d2917SAlp Dener   *accept = PETSC_FALSE;
82128017e9fSAlp Dener   switch (updateType) {
822080d2917SAlp Dener   case BNK_UPDATE_STEP:
823c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
824080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
8259566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step));
826080d2917SAlp Dener       if (step < bnk->nu1) {
827080d2917SAlp Dener         /* Very bad step taken; reduce radius */
828080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
829080d2917SAlp Dener       } else if (step < bnk->nu2) {
830080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
831080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
832080d2917SAlp Dener       } else if (step < bnk->nu3) {
833080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
834080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
835080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
836080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
837080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
838080d2917SAlp Dener         }
839080d2917SAlp Dener       } else if (step < bnk->nu4) {
840080d2917SAlp Dener         /*  Full step taken; increase the radius */
841080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
842080d2917SAlp Dener       } else {
843080d2917SAlp Dener         /*  More than full step taken; increase the radius */
844080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
845080d2917SAlp Dener       }
846080d2917SAlp Dener     } else {
847080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
848080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
849080d2917SAlp Dener     }
850080d2917SAlp Dener     break;
851080d2917SAlp Dener 
852080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
853080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
854e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
855fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
856fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
857fed79b8eSAlp Dener            be rejected! */
858080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8593105154fSTodd Munson       } else {
860b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
861080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
862080d2917SAlp Dener         } else {
8633105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon)) {
864080d2917SAlp Dener             kappa = 1.0;
8653105154fSTodd Munson           } else {
866080d2917SAlp Dener             kappa = actred / prered;
867080d2917SAlp Dener           }
868fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
869080d2917SAlp Dener           if (kappa < bnk->eta1) {
870fed79b8eSAlp Dener             /* Reject the step */
871080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8723105154fSTodd Munson           } else {
873fed79b8eSAlp Dener             /* Accept the step */
874c133c014SAlp Dener             *accept = PETSC_TRUE;
875c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8768d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
877080d2917SAlp Dener               if (kappa < bnk->eta2) {
878080d2917SAlp Dener                 /* Marginal bad step */
879c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
8803105154fSTodd Munson               } else if (kappa < bnk->eta3) {
881fed79b8eSAlp Dener                 /* Reasonable step */
882fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
8833105154fSTodd Munson               } else if (kappa < bnk->eta4) {
884080d2917SAlp Dener                 /* Good step */
885c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
8863105154fSTodd Munson               } else {
887080d2917SAlp Dener                 /* Very good step */
888c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
889080d2917SAlp Dener               }
890c133c014SAlp Dener             }
891080d2917SAlp Dener           }
892080d2917SAlp Dener         }
893080d2917SAlp Dener       }
894080d2917SAlp Dener     } else {
895080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
896080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
897080d2917SAlp Dener     }
898080d2917SAlp Dener     break;
899080d2917SAlp Dener 
900080d2917SAlp Dener   default:
901080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
902b1c2d0e3SAlp Dener       if (prered < 0.0) {
903080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
904080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
905080d2917SAlp Dener         /*  be rejected! */
906080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
907080d2917SAlp Dener       } else {
908b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
909080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
910080d2917SAlp Dener         } else {
911080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
912080d2917SAlp Dener             kappa = 1.0;
913080d2917SAlp Dener           } else {
914080d2917SAlp Dener             kappa = actred / prered;
915080d2917SAlp Dener           }
916080d2917SAlp Dener 
9179566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
918080d2917SAlp Dener           tau_1   = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
919080d2917SAlp Dener           tau_2   = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
920080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
921080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
922080d2917SAlp Dener 
923080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
924080d2917SAlp Dener             /*  Great agreement */
925080d2917SAlp Dener             *accept = PETSC_TRUE;
926080d2917SAlp Dener             if (tau_max < 1.0) {
927080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
928080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
929080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
930080d2917SAlp Dener             } else {
931080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
932080d2917SAlp Dener             }
933080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
934080d2917SAlp Dener             /*  Good agreement */
935080d2917SAlp Dener             *accept = PETSC_TRUE;
936080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
937080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
938080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
939080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
940080d2917SAlp Dener             } else if (tau_max < 1.0) {
941080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
942080d2917SAlp Dener             } else {
943080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
944080d2917SAlp Dener             }
945080d2917SAlp Dener           } else {
946080d2917SAlp Dener             /*  Not good agreement */
947080d2917SAlp Dener             if (tau_min > 1.0) {
948080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
949080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
950080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
951080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
952080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
953080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
954080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
955080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
956080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
957080d2917SAlp Dener             } else {
958080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
959080d2917SAlp Dener             }
960080d2917SAlp Dener           }
961080d2917SAlp Dener         }
962080d2917SAlp Dener       }
963080d2917SAlp Dener     } else {
964080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
965080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
966080d2917SAlp Dener     }
96728017e9fSAlp Dener     break;
968080d2917SAlp Dener   }
969c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
970c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
971fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
9723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
973080d2917SAlp Dener }
974080d2917SAlp Dener 
975eb910715SAlp Dener /* ---------------------------------------------------------- */
976df278d8fSAlp Dener 
977d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
978d71ae5a4SJacob Faibussowitsch {
97962675beeSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
98062675beeSAlp Dener 
98162675beeSAlp Dener   PetscFunctionBegin;
98262675beeSAlp Dener   switch (stepType) {
983d71ae5a4SJacob Faibussowitsch   case BNK_NEWTON:
984d71ae5a4SJacob Faibussowitsch     ++bnk->newt;
985d71ae5a4SJacob Faibussowitsch     break;
986d71ae5a4SJacob Faibussowitsch   case BNK_BFGS:
987d71ae5a4SJacob Faibussowitsch     ++bnk->bfgs;
988d71ae5a4SJacob Faibussowitsch     break;
989d71ae5a4SJacob Faibussowitsch   case BNK_SCALED_GRADIENT:
990d71ae5a4SJacob Faibussowitsch     ++bnk->sgrad;
991d71ae5a4SJacob Faibussowitsch     break;
992d71ae5a4SJacob Faibussowitsch   case BNK_GRADIENT:
993d71ae5a4SJacob Faibussowitsch     ++bnk->grad;
994d71ae5a4SJacob Faibussowitsch     break;
995d71ae5a4SJacob Faibussowitsch   default:
996d71ae5a4SJacob Faibussowitsch     break;
99762675beeSAlp Dener   }
9983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99962675beeSAlp Dener }
100062675beeSAlp Dener 
100162675beeSAlp Dener /* ---------------------------------------------------------- */
100262675beeSAlp Dener 
1003d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetUp_BNK(Tao tao)
1004d71ae5a4SJacob Faibussowitsch {
1005eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
1006e031d6f5SAlp Dener   PetscInt i;
1007eb910715SAlp Dener 
1008eb910715SAlp Dener   PetscFunctionBegin;
100948a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
101048a46eb9SPierre Jolivet   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
101148a46eb9SPierre Jolivet   if (!bnk->W) PetscCall(VecDuplicate(tao->solution, &bnk->W));
101248a46eb9SPierre Jolivet   if (!bnk->Xold) PetscCall(VecDuplicate(tao->solution, &bnk->Xold));
101348a46eb9SPierre Jolivet   if (!bnk->Gold) PetscCall(VecDuplicate(tao->solution, &bnk->Gold));
101448a46eb9SPierre Jolivet   if (!bnk->Xwork) PetscCall(VecDuplicate(tao->solution, &bnk->Xwork));
101548a46eb9SPierre Jolivet   if (!bnk->Gwork) PetscCall(VecDuplicate(tao->solution, &bnk->Gwork));
101648a46eb9SPierre Jolivet   if (!bnk->unprojected_gradient) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient));
101748a46eb9SPierre Jolivet   if (!bnk->unprojected_gradient_old) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient_old));
101848a46eb9SPierre Jolivet   if (!bnk->Diag_min) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_min));
101948a46eb9SPierre Jolivet   if (!bnk->Diag_max) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_max));
1020e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1021c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1022c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
10239566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old)));
10249566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old));
102589da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
10269566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient)));
10279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient));
1028c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
10299566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->Gold)));
10309566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->G_old));
1031c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
10329566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(tao->gradient)));
10339566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg->gradient));
1034c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
10359566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(tao->stepdirection)));
10369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg->stepdirection));
1037c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
10389566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(bnk->bncg, tao->solution));
1039c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
10409566063dSJacob Faibussowitsch     PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its));
10419566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol));
10429566063dSJacob Faibussowitsch     PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin));
10439566063dSJacob Faibussowitsch     PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP));
10449566063dSJacob Faibussowitsch     PetscCall(TaoSetObjective(bnk->bncg, tao->ops->computeobjective, tao->user_objP));
10459566063dSJacob Faibussowitsch     PetscCall(TaoSetGradient(bnk->bncg, NULL, tao->ops->computegradient, tao->user_gradP));
10469566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(bnk->bncg, NULL, tao->ops->computeobjectiveandgradient, tao->user_objgradP));
10479566063dSJacob Faibussowitsch     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg)));
1048c4b75bccSAlp Dener     for (i = 0; i < tao->numbermonitors; ++i) {
10499566063dSJacob Faibussowitsch       PetscCall(TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
10509566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(tao->monitorcontext[i])));
1051e031d6f5SAlp Dener     }
1052e031d6f5SAlp Dener   }
105383c8fe1dSLisandro Dalcin   bnk->X_inactive    = NULL;
105483c8fe1dSLisandro Dalcin   bnk->G_inactive    = NULL;
105583c8fe1dSLisandro Dalcin   bnk->inactive_work = NULL;
105683c8fe1dSLisandro Dalcin   bnk->active_work   = NULL;
105783c8fe1dSLisandro Dalcin   bnk->inactive_idx  = NULL;
105883c8fe1dSLisandro Dalcin   bnk->active_idx    = NULL;
105983c8fe1dSLisandro Dalcin   bnk->active_lower  = NULL;
106083c8fe1dSLisandro Dalcin   bnk->active_upper  = NULL;
106183c8fe1dSLisandro Dalcin   bnk->active_fixed  = NULL;
106283c8fe1dSLisandro Dalcin   bnk->M             = NULL;
106383c8fe1dSLisandro Dalcin   bnk->H_inactive    = NULL;
106483c8fe1dSLisandro Dalcin   bnk->Hpre_inactive = NULL;
10653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1066eb910715SAlp Dener }
1067eb910715SAlp Dener 
1068eb910715SAlp Dener /*------------------------------------------------------------*/
1069df278d8fSAlp Dener 
1070d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDestroy_BNK(Tao tao)
1071d71ae5a4SJacob Faibussowitsch {
1072eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
1073eb910715SAlp Dener 
1074eb910715SAlp Dener   PetscFunctionBegin;
10759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->W));
10769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Xold));
10779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Gold));
10789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Xwork));
10799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Gwork));
10809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->unprojected_gradient));
10819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->unprojected_gradient_old));
10829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Diag_min));
10839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Diag_max));
10849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_lower));
10859566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_upper));
10869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_fixed));
10879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_idx));
10889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->inactive_idx));
10899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
10909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
10919566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&bnk->bncg));
1092a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
10939566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
10943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1095eb910715SAlp Dener }
1096eb910715SAlp Dener 
1097eb910715SAlp Dener /*------------------------------------------------------------*/
1098df278d8fSAlp Dener 
1099d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetFromOptions_BNK(Tao tao, PetscOptionItems *PetscOptionsObject)
1100d71ae5a4SJacob Faibussowitsch {
1101eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
1102eb910715SAlp Dener 
1103eb910715SAlp Dener   PetscFunctionBegin;
1104d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton-Krylov method for bound constrained optimization");
11059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL));
11069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL));
11079566063dSJacob 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));
11089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval, NULL));
11099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin, NULL));
11109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax, NULL));
11119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac, NULL));
11129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin, NULL));
11139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax, NULL));
11149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac, NULL));
11159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac, NULL));
11169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac, NULL));
11179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac, NULL));
11189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1, NULL));
11199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2, NULL));
11209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3, NULL));
11219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4, NULL));
11229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1, NULL));
11239566063dSJacob 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));
11249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3, NULL));
11259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4, NULL));
11269566063dSJacob 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));
11279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1, NULL));
11289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2, NULL));
11299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3, NULL));
11309566063dSJacob 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));
11319566063dSJacob 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));
11329566063dSJacob 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));
11339566063dSJacob 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));
11349566063dSJacob 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));
11359566063dSJacob 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));
11369566063dSJacob 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));
11379566063dSJacob 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));
11389566063dSJacob 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));
11399566063dSJacob 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));
11409566063dSJacob 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));
11419566063dSJacob 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));
11429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i, NULL));
11439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1, NULL));
11449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2, NULL));
11459566063dSJacob 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));
11469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2, NULL));
11479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3, NULL));
11489566063dSJacob 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));
11499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta, NULL));
11509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius, NULL));
11519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius, NULL));
11529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
11539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
11549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
11559566063dSJacob 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));
1156d0609cedSBarry Smith   PetscOptionsHeadEnd();
11578ebe3e4eSStefano Zampini 
11589566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)(tao))->prefix));
11599566063dSJacob Faibussowitsch   PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_cg_"));
11609566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(bnk->bncg));
11618ebe3e4eSStefano Zampini 
11629566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, ((PetscObject)(tao))->prefix));
11639566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_bnk_"));
11649566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
11653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1166eb910715SAlp Dener }
1167eb910715SAlp Dener 
1168eb910715SAlp Dener /*------------------------------------------------------------*/
1169df278d8fSAlp Dener 
1170d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1171d71ae5a4SJacob Faibussowitsch {
1172eb910715SAlp Dener   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
1173eb910715SAlp Dener   PetscInt  nrejects;
1174eb910715SAlp Dener   PetscBool isascii;
1175eb910715SAlp Dener 
1176eb910715SAlp Dener   PetscFunctionBegin;
11779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1178eb910715SAlp Dener   if (isascii) {
11799566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1180b17ffb64SBarry Smith     PetscCall(TaoView(bnk->bncg, viewer));
1181b9ac7092SAlp Dener     if (bnk->M) {
11829566063dSJacob Faibussowitsch       PetscCall(MatLMVMGetRejectCount(bnk->M, &nrejects));
118363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n", nrejects));
1184eb910715SAlp Dener     }
118563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its));
118663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt));
118748a46eb9SPierre Jolivet     if (bnk->M) PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs));
118863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad));
118963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad));
11909566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n"));
119163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  atol: %" PetscInt_FMT "\n", bnk->ksp_atol));
119263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol));
119363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol));
119463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  negc: %" PetscInt_FMT "\n", bnk->ksp_negc));
119563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol));
119663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  iter: %" PetscInt_FMT "\n", bnk->ksp_iter));
119763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  othr: %" PetscInt_FMT "\n", bnk->ksp_othr));
11989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1199eb910715SAlp Dener   }
12003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1201eb910715SAlp Dener }
1202eb910715SAlp Dener 
1203eb910715SAlp Dener /* ---------------------------------------------------------- */
1204df278d8fSAlp Dener 
1205eb910715SAlp Dener /*MC
1206eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
120766ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
12081cb3f120SPierre Jolivet   system of equations to obtain the step direction dk:
1209eb910715SAlp Dener               Hk dk = -gk
12102b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12112b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1212eb910715SAlp Dener 
1213eb910715SAlp Dener     Options Database Keys:
12149fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
12159fa2b5dcSStefano Zampini . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
12169fa2b5dcSStefano Zampini . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
12179fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
12189fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
12199fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
12209fa2b5dcSStefano Zampini . -tao_bnk_sval - (developer) Hessian perturbation starting value
12219fa2b5dcSStefano Zampini . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
12229fa2b5dcSStefano Zampini . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
12239fa2b5dcSStefano Zampini . -tao_bnk_pmin - (developer) minimum Hessian perturbation
12249fa2b5dcSStefano Zampini . -tao_bnk_pmax - (developer) aximum Hessian perturbation
12259fa2b5dcSStefano Zampini . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
12269fa2b5dcSStefano Zampini . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
12279fa2b5dcSStefano Zampini . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
12289fa2b5dcSStefano Zampini . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
12299fa2b5dcSStefano Zampini . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
12309fa2b5dcSStefano Zampini . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction)
12319fa2b5dcSStefano Zampini . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
12329fa2b5dcSStefano Zampini . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
12339fa2b5dcSStefano Zampini . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction)
12349fa2b5dcSStefano Zampini . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
12359fa2b5dcSStefano Zampini . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
12369fa2b5dcSStefano Zampini . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
12379fa2b5dcSStefano Zampini . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
12389fa2b5dcSStefano Zampini . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
12399fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
12409fa2b5dcSStefano Zampini . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
12419fa2b5dcSStefano Zampini . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation)
12429fa2b5dcSStefano Zampini . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
12439fa2b5dcSStefano Zampini . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
12449fa2b5dcSStefano Zampini . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
12459fa2b5dcSStefano Zampini . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
12469fa2b5dcSStefano Zampini . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation)
12479fa2b5dcSStefano Zampini . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step)
12489fa2b5dcSStefano Zampini . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
12499fa2b5dcSStefano Zampini . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step)
12509fa2b5dcSStefano Zampini . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step)
12519fa2b5dcSStefano Zampini . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
12529fa2b5dcSStefano Zampini . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
12539fa2b5dcSStefano Zampini . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step)
12549fa2b5dcSStefano Zampini . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
12559fa2b5dcSStefano Zampini . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
12569fa2b5dcSStefano Zampini . -tao_bnk_mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
12579fa2b5dcSStefano Zampini . -tao_bnk_mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
12589fa2b5dcSStefano Zampini . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
12599fa2b5dcSStefano Zampini . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
12609fa2b5dcSStefano Zampini . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
12619fa2b5dcSStefano Zampini . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
12629fa2b5dcSStefano Zampini - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1263eb910715SAlp Dener 
1264eb910715SAlp Dener   Level: beginner
1265eb910715SAlp Dener M*/
1266eb910715SAlp Dener 
1267d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoCreate_BNK(Tao tao)
1268d71ae5a4SJacob Faibussowitsch {
1269eb910715SAlp Dener   TAO_BNK *bnk;
1270b9ac7092SAlp Dener   PC       pc;
1271eb910715SAlp Dener 
1272eb910715SAlp Dener   PetscFunctionBegin;
12734dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&bnk));
1274eb910715SAlp Dener 
1275eb910715SAlp Dener   tao->ops->setup          = TaoSetUp_BNK;
1276eb910715SAlp Dener   tao->ops->view           = TaoView_BNK;
1277eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1278eb910715SAlp Dener   tao->ops->destroy        = TaoDestroy_BNK;
1279eb910715SAlp Dener 
1280eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1281eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1282eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1283eb910715SAlp Dener 
1284eb910715SAlp Dener   tao->data = (void *)bnk;
1285eb910715SAlp Dener 
128666ed3702SAlp Dener   /*  Hessian shifting parameters */
1287e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1288e0ed867bSAlp Dener   bnk->computestep    = TaoBNKComputeStep;
1289e0ed867bSAlp Dener 
1290eb910715SAlp Dener   bnk->sval  = 0.0;
1291eb910715SAlp Dener   bnk->imin  = 1.0e-4;
1292eb910715SAlp Dener   bnk->imax  = 1.0e+2;
1293eb910715SAlp Dener   bnk->imfac = 1.0e-1;
1294eb910715SAlp Dener 
1295eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1296eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1297eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1298eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1299eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1300eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1301eb910715SAlp Dener 
1302eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1303eb910715SAlp Dener   bnk->nu1 = 0.25;
1304eb910715SAlp Dener   bnk->nu2 = 0.50;
1305eb910715SAlp Dener   bnk->nu3 = 1.00;
1306eb910715SAlp Dener   bnk->nu4 = 1.25;
1307eb910715SAlp Dener 
1308eb910715SAlp Dener   bnk->omega1 = 0.25;
1309eb910715SAlp Dener   bnk->omega2 = 0.50;
1310eb910715SAlp Dener   bnk->omega3 = 1.00;
1311eb910715SAlp Dener   bnk->omega4 = 2.00;
1312eb910715SAlp Dener   bnk->omega5 = 4.00;
1313eb910715SAlp Dener 
1314eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1315eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1316eb910715SAlp Dener   bnk->eta2 = 0.25;
1317eb910715SAlp Dener   bnk->eta3 = 0.50;
1318eb910715SAlp Dener   bnk->eta4 = 0.90;
1319eb910715SAlp Dener 
1320eb910715SAlp Dener   bnk->alpha1 = 0.25;
1321eb910715SAlp Dener   bnk->alpha2 = 0.50;
1322eb910715SAlp Dener   bnk->alpha3 = 1.00;
1323eb910715SAlp Dener   bnk->alpha4 = 2.00;
1324eb910715SAlp Dener   bnk->alpha5 = 4.00;
1325eb910715SAlp Dener 
1326eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1327eb910715SAlp Dener   bnk->mu1 = 0.10;
1328eb910715SAlp Dener   bnk->mu2 = 0.50;
1329eb910715SAlp Dener 
1330eb910715SAlp Dener   bnk->gamma1 = 0.25;
1331eb910715SAlp Dener   bnk->gamma2 = 0.50;
1332eb910715SAlp Dener   bnk->gamma3 = 2.00;
1333eb910715SAlp Dener   bnk->gamma4 = 4.00;
1334eb910715SAlp Dener 
1335eb910715SAlp Dener   bnk->theta = 0.05;
1336eb910715SAlp Dener 
1337eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1338eb910715SAlp Dener   bnk->mu1_i = 0.35;
1339eb910715SAlp Dener   bnk->mu2_i = 0.50;
1340eb910715SAlp Dener 
1341eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1342eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1343eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1344eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1345eb910715SAlp Dener 
1346eb910715SAlp Dener   bnk->theta_i = 0.25;
1347eb910715SAlp Dener 
1348eb910715SAlp Dener   /*  Remaining parameters */
1349c0f10754SAlp Dener   bnk->max_cg_its = 0;
1350eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1351eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1352770b7498SAlp Dener   bnk->epsilon    = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0);
13530a4511e9SAlp Dener   bnk->as_tol     = 1.0e-3;
13540a4511e9SAlp Dener   bnk->as_step    = 1.0e-3;
135562675beeSAlp Dener   bnk->dmin       = 1.0e-6;
135662675beeSAlp Dener   bnk->dmax       = 1.0e6;
1357eb910715SAlp Dener 
135883c8fe1dSLisandro Dalcin   bnk->M           = NULL;
135983c8fe1dSLisandro Dalcin   bnk->bfgs_pre    = NULL;
1360eb910715SAlp Dener   bnk->init_type   = BNK_INIT_INTERPOLATION;
13617b1c7716SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION;
13622f75a4aaSAlp Dener   bnk->as_type     = BNK_AS_BERTSEKAS;
1363eb910715SAlp Dener 
1364e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
13659566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg));
13669566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1));
13679566063dSJacob Faibussowitsch   PetscCall(TaoSetType(bnk->bncg, TAOBNCG));
1368e031d6f5SAlp Dener 
1369c0f10754SAlp Dener   /* Create the line search */
13709566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
13719566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1372f4db9bf7SStefano Zampini   PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
13739566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
1374eb910715SAlp Dener 
1375eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
13769566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
13779566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
13789566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
13799566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
13809566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLMVM));
13813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1382eb910715SAlp Dener }
1383