xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener #include <petscksp.h>
4eb910715SAlp Dener 
570a3f44bSAlp Dener static const char *BNK_INIT[64]   = {"constant", "direction", "interpolation"};
670a3f44bSAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
770a3f44bSAlp Dener static const char *BNK_AS[64]     = {"none", "bertsekas"};
870a3f44bSAlp Dener 
9e031d6f5SAlp Dener /*------------------------------------------------------------*/
10e031d6f5SAlp Dener 
11df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
12df278d8fSAlp Dener 
139371c9d4SSatish Balay PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH) {
14eb910715SAlp Dener   TAO_BNK          *bnk = (TAO_BNK *)tao->data;
15eb910715SAlp Dener   PC                pc;
1689da521bSAlp Dener   PetscReal         f_min, ftrial, prered, actred, kappa, sigma, resnorm;
17eb910715SAlp Dener   PetscReal         tau, tau_1, tau_2, tau_max, tau_min, max_radius;
180ad3a497SAlp Dener   PetscBool         is_bfgs, is_jacobi, is_symmetric, sym_set;
19c4b75bccSAlp Dener   PetscInt          n, N, nDiff;
20eb910715SAlp Dener   PetscInt          i_max = 5;
21eb910715SAlp Dener   PetscInt          j_max = 1;
22eb910715SAlp Dener   PetscInt          i, j;
232e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
24eb910715SAlp Dener 
25eb910715SAlp Dener   PetscFunctionBegin;
2628017e9fSAlp Dener   /* Project the current point onto the feasible set */
279566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
289566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU));
291baa6e33SBarry Smith   if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
3028017e9fSAlp Dener 
3128017e9fSAlp Dener   /* Project the initial point onto the feasible region */
329566063dSJacob Faibussowitsch   PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
3328017e9fSAlp Dener 
3428017e9fSAlp Dener   /* Check convergence criteria */
359566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
369566063dSJacob Faibussowitsch   PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
379566063dSJacob Faibussowitsch   PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
389566063dSJacob Faibussowitsch   PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
399566063dSJacob Faibussowitsch   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
4028017e9fSAlp Dener 
41c0f10754SAlp Dener   /* Test the initial point for convergence */
429566063dSJacob Faibussowitsch   PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
439566063dSJacob Faibussowitsch   PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
443c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
459566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
469566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
47dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
48c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
49c0f10754SAlp Dener 
50e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
51eb910715SAlp Dener   bnk->ksp_atol = 0;
52eb910715SAlp Dener   bnk->ksp_rtol = 0;
53eb910715SAlp Dener   bnk->ksp_dtol = 0;
54eb910715SAlp Dener   bnk->ksp_ctol = 0;
55eb910715SAlp Dener   bnk->ksp_negc = 0;
56eb910715SAlp Dener   bnk->ksp_iter = 0;
57eb910715SAlp Dener   bnk->ksp_othr = 0;
58eb910715SAlp Dener 
59e031d6f5SAlp Dener   /* Reset accepted step type counters */
60e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
61e031d6f5SAlp Dener   bnk->newt       = 0;
62e031d6f5SAlp Dener   bnk->bfgs       = 0;
63e031d6f5SAlp Dener   bnk->sgrad      = 0;
64e031d6f5SAlp Dener   bnk->grad       = 0;
65e031d6f5SAlp Dener 
66fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
67fed79b8eSAlp Dener   bnk->pert = bnk->sval;
68fed79b8eSAlp Dener 
69937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
709566063dSJacob Faibussowitsch   PetscCall(VecSet(tao->stepdirection, 0.0));
71937a31a1SAlp Dener 
72e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
739566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
76b9ac7092SAlp Dener   if (is_bfgs) {
77b9ac7092SAlp Dener     bnk->bfgs_pre = pc;
789566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M));
799566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->solution, &n));
809566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->solution, &N));
819566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(bnk->M, n, n, N, N));
829566063dSJacob Faibussowitsch     PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient));
839566063dSJacob Faibussowitsch     PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric));
843c859ba3SBarry Smith     PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
851baa6e33SBarry Smith   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
86e031d6f5SAlp Dener 
87e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
889566063dSJacob Faibussowitsch   PetscCall(VecSet(bnk->Diag_min, bnk->dmin));
899566063dSJacob Faibussowitsch   PetscCall(VecSet(bnk->Diag_max, bnk->dmax));
90eb910715SAlp Dener 
91eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
92eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
93c0f10754SAlp Dener   *needH = PETSC_TRUE;
949566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp, "KSPCGSetRadius_C", &kspTR));
952e6e4ca1SStefano Zampini   if (kspTR) {
9662675beeSAlp Dener     switch (initType) {
97eb910715SAlp Dener     case BNK_INIT_CONSTANT:
98eb910715SAlp Dener       /* Use the initial radius specified */
99c0f10754SAlp Dener       tao->trust = tao->trust0;
100eb910715SAlp Dener       break;
101eb910715SAlp Dener 
102eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
103c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
104eb910715SAlp Dener       max_radius = 0.0;
10508752603SAlp Dener       tao->trust = tao->trust0;
106eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1070a4511e9SAlp Dener         f_min = bnk->f;
108eb910715SAlp Dener         sigma = 0.0;
109eb910715SAlp Dener 
110c0f10754SAlp Dener         if (*needH) {
11162602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
1129566063dSJacob Faibussowitsch           PetscCall((*bnk->computehessian)(tao));
1139566063dSJacob Faibussowitsch           PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE));
1149566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&bnk->H_inactive));
11589da521bSAlp Dener           if (bnk->active_idx) {
1169566063dSJacob Faibussowitsch             PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
11728017e9fSAlp Dener           } else {
1189566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)tao->hessian));
119c5e9d94cSAlp Dener             bnk->H_inactive = tao->hessian;
12028017e9fSAlp Dener           }
121c0f10754SAlp Dener           *needH = PETSC_FALSE;
122eb910715SAlp Dener         }
123eb910715SAlp Dener 
124eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
12562602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
1269566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->Xold));
1279566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, -tao->trust / bnk->gnorm, tao->gradient));
1289566063dSJacob Faibussowitsch           PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
12989da521bSAlp Dener           /* Compute the step we actually accepted */
1309566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->W));
1319566063dSJacob Faibussowitsch           PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold));
13262602cfbSAlp Dener           /* Compute the objective at the trial */
1339566063dSJacob Faibussowitsch           PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial));
1343c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
1359566063dSJacob Faibussowitsch           PetscCall(VecCopy(bnk->Xold, tao->solution));
136eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
137eb910715SAlp Dener             tau = bnk->gamma1_i;
138eb910715SAlp Dener           } else {
1390a4511e9SAlp Dener             if (ftrial < f_min) {
1400a4511e9SAlp Dener               f_min = ftrial;
141eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
142eb910715SAlp Dener             }
14308752603SAlp Dener 
144770b7498SAlp Dener             /* Compute the predicted and actual reduction */
14589da521bSAlp Dener             if (bnk->active_idx) {
1469566063dSJacob Faibussowitsch               PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1479566063dSJacob Faibussowitsch               PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1482ab2a32cSAlp Dener             } else {
14908752603SAlp Dener               bnk->X_inactive    = bnk->W;
15008752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1512ab2a32cSAlp Dener             }
1529566063dSJacob Faibussowitsch             PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
1539566063dSJacob Faibussowitsch             PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered));
15489da521bSAlp Dener             if (bnk->active_idx) {
1559566063dSJacob Faibussowitsch               PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1569566063dSJacob Faibussowitsch               PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1572ab2a32cSAlp Dener             }
158eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
159eb910715SAlp Dener             actred = bnk->f - ftrial;
1603105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
161eb910715SAlp Dener               kappa = 1.0;
1623105154fSTodd Munson             } else {
163eb910715SAlp Dener               kappa = actred / prered;
164eb910715SAlp Dener             }
165eb910715SAlp Dener 
166eb910715SAlp Dener             tau_1   = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
167eb910715SAlp Dener             tau_2   = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
168eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
169eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
170eb910715SAlp Dener 
17118cfbf8eSSatish Balay             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
172eb910715SAlp Dener               /*  Great agreement */
173eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
174eb910715SAlp Dener 
175eb910715SAlp Dener               if (tau_max < 1.0) {
176eb910715SAlp Dener                 tau = bnk->gamma3_i;
1773105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
178eb910715SAlp Dener                 tau = bnk->gamma4_i;
1793105154fSTodd Munson               } else {
180eb910715SAlp Dener                 tau = tau_max;
181eb910715SAlp Dener               }
18218cfbf8eSSatish Balay             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
183eb910715SAlp Dener               /*  Good agreement */
184eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
185eb910715SAlp Dener 
186eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
187eb910715SAlp Dener                 tau = bnk->gamma2_i;
188eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
189eb910715SAlp Dener                 tau = bnk->gamma3_i;
190eb910715SAlp Dener               } else {
191eb910715SAlp Dener                 tau = tau_max;
192eb910715SAlp Dener               }
1938f8a4e06SAlp Dener             } else {
194eb910715SAlp Dener               /*  Not good agreement */
195eb910715SAlp Dener               if (tau_min > 1.0) {
196eb910715SAlp Dener                 tau = bnk->gamma2_i;
197eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
198eb910715SAlp Dener                 tau = bnk->gamma1_i;
199eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
200eb910715SAlp Dener                 tau = bnk->gamma1_i;
2013105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
202eb910715SAlp Dener                 tau = tau_1;
2033105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
204eb910715SAlp Dener                 tau = tau_2;
205eb910715SAlp Dener               } else {
206eb910715SAlp Dener                 tau = tau_max;
207eb910715SAlp Dener               }
208eb910715SAlp Dener             }
209eb910715SAlp Dener           }
210eb910715SAlp Dener           tao->trust = tau * tao->trust;
211eb910715SAlp Dener         }
212eb910715SAlp Dener 
2130a4511e9SAlp Dener         if (f_min < bnk->f) {
214937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2150a4511e9SAlp Dener           bnk->f = f_min;
2169566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, bnk->Xold));
2179566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
2189566063dSJacob Faibussowitsch           PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
2199566063dSJacob Faibussowitsch           PetscCall(VecCopy(tao->solution, tao->stepdirection));
2209566063dSJacob Faibussowitsch           PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
2219566063dSJacob Faibussowitsch           PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
2229566063dSJacob Faibussowitsch           PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
2239566063dSJacob Faibussowitsch           PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
2249566063dSJacob Faibussowitsch           PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
225937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
2269566063dSJacob Faibussowitsch           PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
227c0f10754SAlp Dener           *needH = PETSC_TRUE;
228937a31a1SAlp Dener           /* Test the new step for convergence */
2299566063dSJacob Faibussowitsch           PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
2309566063dSJacob Faibussowitsch           PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
2313c859ba3SBarry Smith           PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2329566063dSJacob Faibussowitsch           PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
2339566063dSJacob Faibussowitsch           PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
234dbbe0bcdSBarry Smith           PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
235eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
236937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
2379566063dSJacob Faibussowitsch           PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
238eb910715SAlp Dener         }
239eb910715SAlp Dener       }
240eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
241e031d6f5SAlp Dener 
242e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
243e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
244e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
245eb910715SAlp Dener       break;
246eb910715SAlp Dener 
247eb910715SAlp Dener     default:
248eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
249eb910715SAlp Dener       tao->trust = 0.0;
250eb910715SAlp Dener       break;
251eb910715SAlp Dener     }
252eb910715SAlp Dener   }
253eb910715SAlp Dener   PetscFunctionReturn(0);
254eb910715SAlp Dener }
255eb910715SAlp Dener 
256df278d8fSAlp Dener /*------------------------------------------------------------*/
257df278d8fSAlp Dener 
258e0ed867bSAlp Dener /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */
25962675beeSAlp Dener 
2609371c9d4SSatish Balay PetscErrorCode TaoBNKComputeHessian(Tao tao) {
26162675beeSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
26262675beeSAlp Dener 
26362675beeSAlp Dener   PetscFunctionBegin;
26462675beeSAlp Dener   /* Compute the Hessian */
2659566063dSJacob Faibussowitsch   PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
26662675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
2671baa6e33SBarry Smith   if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
268e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
2699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
2709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
271f5766c09SAlp Dener   if (bnk->active_idx) {
2729566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
273e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
2749566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
275e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
276e0ed867bSAlp Dener     } else {
2779566063dSJacob Faibussowitsch       PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive));
278e0ed867bSAlp Dener     }
2791baa6e33SBarry Smith     if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx));
280e0ed867bSAlp Dener   } else {
2819566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)tao->hessian));
282c5e9d94cSAlp Dener     bnk->H_inactive = tao->hessian;
283e0ed867bSAlp Dener     if (tao->hessian == tao->hessian_pre) {
2849566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
285e0ed867bSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
286e0ed867bSAlp Dener     } else {
2879566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre));
288c5e9d94cSAlp Dener       bnk->Hpre_inactive = tao->hessian_pre;
289e0ed867bSAlp Dener     }
2901baa6e33SBarry Smith     if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre));
291e0ed867bSAlp Dener   }
29262675beeSAlp Dener   PetscFunctionReturn(0);
29362675beeSAlp Dener }
29462675beeSAlp Dener 
29562675beeSAlp Dener /*------------------------------------------------------------*/
29662675beeSAlp Dener 
2972f75a4aaSAlp Dener /* Routine for estimating the active set */
2982f75a4aaSAlp Dener 
2999371c9d4SSatish Balay PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType) {
3002f75a4aaSAlp Dener   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
301f4db9bf7SStefano Zampini   PetscBool hessComputed, diagExists, hadactive;
3022f75a4aaSAlp Dener 
3032f75a4aaSAlp Dener   PetscFunctionBegin;
304f4db9bf7SStefano Zampini   hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE;
30508752603SAlp Dener   switch (asType) {
3062f75a4aaSAlp Dener   case BNK_AS_NONE:
3079566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bnk->inactive_idx));
3089566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx));
3099566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bnk->active_idx));
3109566063dSJacob Faibussowitsch     PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx));
3112f75a4aaSAlp Dener     break;
3122f75a4aaSAlp Dener 
3132f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3142f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
315b9ac7092SAlp Dener     if (bnk->M) {
3162f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3179566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W));
3182f75a4aaSAlp Dener     } else {
319fc5ca067SStefano Zampini       hessComputed = diagExists = PETSC_FALSE;
320*48a46eb9SPierre Jolivet       if (tao->hessian) PetscCall(MatAssembled(tao->hessian, &hessComputed));
321*48a46eb9SPierre Jolivet       if (hessComputed) PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists));
322fc5ca067SStefano Zampini       if (diagExists) {
3239b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3249566063dSJacob Faibussowitsch         PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork));
3259566063dSJacob Faibussowitsch         PetscCall(VecAbs(bnk->Xwork));
3269566063dSJacob Faibussowitsch         PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork));
3279566063dSJacob Faibussowitsch         PetscCall(VecReciprocal(bnk->Xwork));
3289566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient));
32961be54a6SAlp Dener       } else {
330c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
3319566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W));
33261be54a6SAlp Dener       }
3332f75a4aaSAlp Dener     }
3349566063dSJacob Faibussowitsch     PetscCall(VecScale(bnk->W, -1.0));
3359371c9d4SSatish 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));
336c4b75bccSAlp Dener     break;
3372f75a4aaSAlp Dener 
3389371c9d4SSatish Balay   default: break;
3392f75a4aaSAlp Dener   }
340f4db9bf7SStefano Zampini   bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */
3412f75a4aaSAlp Dener   PetscFunctionReturn(0);
3422f75a4aaSAlp Dener }
3432f75a4aaSAlp Dener 
3442f75a4aaSAlp Dener /*------------------------------------------------------------*/
3452f75a4aaSAlp Dener 
3462f75a4aaSAlp Dener /* Routine for bounding the step direction */
3472f75a4aaSAlp Dener 
3489371c9d4SSatish Balay PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step) {
3492f75a4aaSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
3502f75a4aaSAlp Dener 
3512f75a4aaSAlp Dener   PetscFunctionBegin;
352a1318120SAlp Dener   switch (asType) {
3539371c9d4SSatish Balay   case BNK_AS_NONE: PetscCall(VecISSet(step, bnk->active_idx, 0.0)); break;
3542f75a4aaSAlp Dener 
3559371c9d4SSatish Balay   case BNK_AS_BERTSEKAS: PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step)); break;
3562f75a4aaSAlp Dener 
3579371c9d4SSatish Balay   default: break;
3582f75a4aaSAlp Dener   }
3592f75a4aaSAlp Dener   PetscFunctionReturn(0);
3602f75a4aaSAlp Dener }
3612f75a4aaSAlp Dener 
362e031d6f5SAlp Dener /*------------------------------------------------------------*/
363e031d6f5SAlp Dener 
364e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
365e031d6f5SAlp Dener    accelerate Newton convergence.
366e031d6f5SAlp Dener 
367e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
368e031d6f5SAlp Dener    for more gradient evaluations.
369e031d6f5SAlp Dener */
370e031d6f5SAlp Dener 
3719371c9d4SSatish Balay PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate) {
372c0f10754SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
373c0f10754SAlp Dener 
374c0f10754SAlp Dener   PetscFunctionBegin;
375c0f10754SAlp Dener   *terminate = PETSC_FALSE;
376c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
377c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
378c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
379c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
3809566063dSJacob Faibussowitsch     PetscCall(TaoSolve(bnk->bncg));
381c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
382c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
383c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
384c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
385c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
386e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
387c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
388c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
389c0f10754SAlp 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) {
390c0f10754SAlp Dener       *terminate = PETSC_TRUE;
39161be54a6SAlp Dener     } else {
3929566063dSJacob Faibussowitsch       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
393c0f10754SAlp Dener     }
394c0f10754SAlp Dener   }
395c0f10754SAlp Dener   PetscFunctionReturn(0);
396c0f10754SAlp Dener }
397c0f10754SAlp Dener 
3982f75a4aaSAlp Dener /*------------------------------------------------------------*/
3992f75a4aaSAlp Dener 
400c0f10754SAlp Dener /* Routine for computing the Newton step. */
401df278d8fSAlp Dener 
4029371c9d4SSatish Balay PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) {
403eb910715SAlp Dener   TAO_BNK          *bnk         = (TAO_BNK *)tao->data;
404eb910715SAlp Dener   PetscInt          bfgsUpdates = 0;
405eb910715SAlp Dener   PetscInt          kspits;
406bddd1ffdSAlp Dener   PetscBool         is_lmvm;
4072e6e4ca1SStefano Zampini   PetscVoidFunction kspTR;
408eb910715SAlp Dener 
409eb910715SAlp Dener   PetscFunctionBegin;
41089da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
41189da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
41289da521bSAlp Dener   if (!bnk->inactive_idx) {
4139566063dSJacob Faibussowitsch     PetscCall(VecSet(tao->stepdirection, 0.0));
4149566063dSJacob Faibussowitsch     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
41589da521bSAlp Dener     PetscFunctionReturn(0);
41689da521bSAlp Dener   }
41789da521bSAlp Dener 
41862675beeSAlp Dener   /* Shift the reduced Hessian matrix */
419e831869dSStefano Zampini   if (shift && bnk->pert > 0) {
4209566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm));
421f7bf01afSAlp Dener     if (is_lmvm) {
4229566063dSJacob Faibussowitsch       PetscCall(MatShift(tao->hessian, bnk->pert));
423f7bf01afSAlp Dener     } else {
4249566063dSJacob Faibussowitsch       PetscCall(MatShift(bnk->H_inactive, bnk->pert));
425*48a46eb9SPierre Jolivet       if (bnk->H_inactive != bnk->Hpre_inactive) PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert));
42662675beeSAlp Dener     }
427f7bf01afSAlp Dener   }
42862675beeSAlp Dener 
429eb910715SAlp Dener   /* Solve the Newton system of equations */
430937a31a1SAlp Dener   tao->ksp_its = 0;
4319566063dSJacob Faibussowitsch   PetscCall(VecSet(tao->stepdirection, 0.0));
432f4db9bf7SStefano Zampini   if (bnk->resetksp) {
4339566063dSJacob Faibussowitsch     PetscCall(KSPReset(tao->ksp));
4349566063dSJacob Faibussowitsch     PetscCall(KSPResetFromOptions(tao->ksp));
435f4db9bf7SStefano Zampini     bnk->resetksp = PETSC_FALSE;
436f4db9bf7SStefano Zampini   }
4379566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(tao->ksp, bnk->H_inactive, bnk->Hpre_inactive));
4389566063dSJacob Faibussowitsch   PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork));
43989da521bSAlp Dener   if (bnk->active_idx) {
4409566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4419566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4425e9b73cbSAlp Dener   } else {
4435e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
4445e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
44528017e9fSAlp Dener   }
4469566063dSJacob Faibussowitsch   PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
4479566063dSJacob Faibussowitsch   PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4489566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
449eb910715SAlp Dener   tao->ksp_its += kspits;
450eb910715SAlp Dener   tao->ksp_tot_its += kspits;
451f4db9bf7SStefano Zampini   PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp, "KSPCGGetNormD_C", &kspTR));
452f4db9bf7SStefano Zampini   if (kspTR) {
4539566063dSJacob Faibussowitsch     PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
454eb910715SAlp Dener 
455eb910715SAlp Dener     if (0.0 == tao->trust) {
456eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
457080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
458080d2917SAlp Dener         tao->trust = bnk->dnorm;
459eb910715SAlp Dener 
460eb910715SAlp Dener         /* Modify the radius if it is too large or small */
461eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
462eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
463eb910715SAlp Dener       } else {
464eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
465eb910715SAlp Dener            the trust-region subproblem to get a direction */
466eb910715SAlp Dener         tao->trust = tao->trust0;
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 
4729566063dSJacob Faibussowitsch         PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
4739566063dSJacob Faibussowitsch         PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4749566063dSJacob Faibussowitsch         PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
475eb910715SAlp Dener         tao->ksp_its += kspits;
476eb910715SAlp Dener         tao->ksp_tot_its += kspits;
4779566063dSJacob Faibussowitsch         PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
478eb910715SAlp Dener 
4793c859ba3SBarry Smith         PetscCheck(bnk->dnorm != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
480eb910715SAlp Dener       }
481eb910715SAlp Dener     }
482eb910715SAlp Dener   }
4835e9b73cbSAlp Dener   /* Restore sub vectors back */
48489da521bSAlp Dener   if (bnk->active_idx) {
4859566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4869566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4875e9b73cbSAlp Dener   }
488770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
4899566063dSJacob Faibussowitsch   PetscCall(VecScale(tao->stepdirection, -1.0));
4909566063dSJacob Faibussowitsch   PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
491770b7498SAlp Dener 
492770b7498SAlp Dener   /* Record convergence reasons */
4939566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason));
494e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
495770b7498SAlp Dener     ++bnk->ksp_atol;
496e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
497770b7498SAlp Dener     ++bnk->ksp_rtol;
498e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
499770b7498SAlp Dener     ++bnk->ksp_ctol;
500e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
501770b7498SAlp Dener     ++bnk->ksp_negc;
502e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
503770b7498SAlp Dener     ++bnk->ksp_dtol;
504e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
505770b7498SAlp Dener     ++bnk->ksp_iter;
506770b7498SAlp Dener   } else {
507770b7498SAlp Dener     ++bnk->ksp_othr;
508770b7498SAlp Dener   }
509fed79b8eSAlp Dener 
510fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
511b9ac7092SAlp Dener   if (bnk->M) {
5129566063dSJacob Faibussowitsch     PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
513b2d8c577SAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
514fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
5159566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
5169566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
517eb910715SAlp Dener     }
518fed79b8eSAlp Dener   }
5196b591159SAlp Dener   *step_type = BNK_NEWTON;
520e465cd6fSAlp Dener   PetscFunctionReturn(0);
521e465cd6fSAlp Dener }
522eb910715SAlp Dener 
52362675beeSAlp Dener /*------------------------------------------------------------*/
52462675beeSAlp Dener 
5255e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5265e9b73cbSAlp Dener 
5279371c9d4SSatish Balay PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered) {
5285e9b73cbSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
5295e9b73cbSAlp Dener 
5305e9b73cbSAlp Dener   PetscFunctionBegin;
5315e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
53289da521bSAlp Dener   if (bnk->active_idx) {
5339566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5349566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5359566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5365e9b73cbSAlp Dener   } else {
5375e9b73cbSAlp Dener     bnk->X_inactive    = tao->stepdirection;
5385e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
5395e9b73cbSAlp Dener     bnk->G_inactive    = bnk->Gwork;
5405e9b73cbSAlp Dener   }
5415e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
5429566063dSJacob Faibussowitsch   PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
5439566063dSJacob Faibussowitsch   PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive));
5449566063dSJacob Faibussowitsch   PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered));
5455e9b73cbSAlp Dener   /* Restore the sub vectors */
54689da521bSAlp Dener   if (bnk->active_idx) {
5479566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5489566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5499566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5505e9b73cbSAlp Dener   }
5515e9b73cbSAlp Dener   PetscFunctionReturn(0);
5525e9b73cbSAlp Dener }
5535e9b73cbSAlp Dener 
5545e9b73cbSAlp Dener /*------------------------------------------------------------*/
5555e9b73cbSAlp Dener 
55662675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
55762675beeSAlp Dener 
55862675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
55962675beeSAlp Dener    in the event that the Newton step fails the test.
56062675beeSAlp Dener */
56162675beeSAlp Dener 
5629371c9d4SSatish Balay PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType) {
563e465cd6fSAlp Dener   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
564b2d8c577SAlp Dener   PetscReal gdx, e_min;
565e465cd6fSAlp Dener   PetscInt  bfgsUpdates;
566e465cd6fSAlp Dener 
567e465cd6fSAlp Dener   PetscFunctionBegin;
5686b591159SAlp Dener   switch (*stepType) {
5696b591159SAlp Dener   case BNK_NEWTON:
5709566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
571eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
572eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
573eb910715SAlp Dener         Update the perturbation for next time */
574eb910715SAlp Dener       if (bnk->pert <= 0.0) {
5752e6e4ca1SStefano Zampini         PetscBool is_gltr;
5762e6e4ca1SStefano Zampini 
577eb910715SAlp Dener         /* Initialize the perturbation */
578eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
5799566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr));
5802e6e4ca1SStefano Zampini         if (is_gltr) {
5819566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
582eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
583eb910715SAlp Dener         }
584eb910715SAlp Dener       } else {
585eb910715SAlp Dener         /* Increase the perturbation */
586eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
587eb910715SAlp Dener       }
588eb910715SAlp Dener 
5890ad3a497SAlp Dener       if (!bnk->M) {
590eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
591eb910715SAlp Dener           Must use gradient direction in this case */
5929566063dSJacob Faibussowitsch         PetscCall(VecCopy(tao->gradient, tao->stepdirection));
593eb910715SAlp Dener         *stepType = BNK_GRADIENT;
594eb910715SAlp Dener       } else {
595eb910715SAlp Dener         /* Attempt to use the BFGS direction */
5969566063dSJacob Faibussowitsch         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
597eb910715SAlp Dener 
5988d5ead36SAlp Dener         /* Check for success (descent direction)
5998d5ead36SAlp Dener           NOTE: Negative gdx here means not a descent direction because
6008d5ead36SAlp Dener           the fall-back step is missing a negative sign. */
6019566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
6023105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
603eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
604eb910715SAlp Dener             We can assert bfgsUpdates > 1 in this case because
605eb910715SAlp Dener             the first solve produces the scaled gradient direction,
606eb910715SAlp Dener             which is guaranteed to be descent */
607eb910715SAlp Dener 
608eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
6099566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6109566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6119566063dSJacob Faibussowitsch           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
612eb910715SAlp Dener 
613eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
614eb910715SAlp Dener         } else {
6159566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
616eb910715SAlp Dener           if (1 == bfgsUpdates) {
617eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
618eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
619eb910715SAlp Dener           } else {
620eb910715SAlp Dener             *stepType = BNK_BFGS;
621eb910715SAlp Dener           }
622eb910715SAlp Dener         }
623eb910715SAlp Dener       }
6248d5ead36SAlp Dener       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6259566063dSJacob Faibussowitsch       PetscCall(VecScale(tao->stepdirection, -1.0));
6269566063dSJacob Faibussowitsch       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
627eb910715SAlp Dener     } else {
628eb910715SAlp Dener       /* Computed Newton step is descent */
629eb910715SAlp Dener       switch (ksp_reason) {
630eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
631eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
632eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
633eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
634eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
635eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
636eb910715SAlp Dener         if (bnk->pert <= 0.0) {
6372e6e4ca1SStefano Zampini           PetscBool is_gltr;
6382e6e4ca1SStefano Zampini 
639eb910715SAlp Dener           /* Initialize the perturbation */
640eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
6419566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr));
6422e6e4ca1SStefano Zampini           if (is_gltr) {
6439566063dSJacob Faibussowitsch             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
644eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
645eb910715SAlp Dener           }
646eb910715SAlp Dener         } else {
647eb910715SAlp Dener           /* Increase the perturbation */
648eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
649eb910715SAlp Dener         }
650eb910715SAlp Dener         break;
651eb910715SAlp Dener 
652eb910715SAlp Dener       default:
653eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
654eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
6559371c9d4SSatish Balay         if (bnk->pert < bnk->pmin) { bnk->pert = 0.0; }
656eb910715SAlp Dener         break;
657eb910715SAlp Dener       }
658fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
659eb910715SAlp Dener     }
6606b591159SAlp Dener     break;
6616b591159SAlp Dener 
6626b591159SAlp Dener   case BNK_BFGS:
6636b591159SAlp Dener     /* Check for success (descent direction) */
6649566063dSJacob Faibussowitsch     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
6656b591159SAlp Dener     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
6666b591159SAlp Dener       /* Step is not descent or solve was not successful
6676b591159SAlp Dener          Use steepest descent direction (scaled) */
6689566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6699566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6709566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection));
6719566063dSJacob Faibussowitsch       PetscCall(VecScale(tao->stepdirection, -1.0));
6729566063dSJacob Faibussowitsch       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
6736b591159SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
6746b591159SAlp Dener     } else {
6756b591159SAlp Dener       *stepType = BNK_BFGS;
6766b591159SAlp Dener     }
6776b591159SAlp Dener     break;
6786b591159SAlp Dener 
6799371c9d4SSatish Balay   case BNK_SCALED_GRADIENT: break;
6806b591159SAlp Dener 
6819371c9d4SSatish Balay   default: break;
6826b591159SAlp Dener   }
6836b591159SAlp Dener 
684eb910715SAlp Dener   PetscFunctionReturn(0);
685eb910715SAlp Dener }
686eb910715SAlp Dener 
687df278d8fSAlp Dener /*------------------------------------------------------------*/
688df278d8fSAlp Dener 
689df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
690df278d8fSAlp Dener 
691df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
692df278d8fSAlp Dener   Newton step does not produce a valid step length.
693df278d8fSAlp Dener */
694df278d8fSAlp Dener 
6959371c9d4SSatish Balay PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) {
696c14b763aSAlp Dener   TAO_BNK                     *bnk = (TAO_BNK *)tao->data;
697c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
698b2d8c577SAlp Dener   PetscReal                    e_min, gdx;
699c14b763aSAlp Dener   PetscInt                     bfgsUpdates;
700c14b763aSAlp Dener 
701c14b763aSAlp Dener   PetscFunctionBegin;
702c14b763aSAlp Dener   /* Perform the linesearch */
7039566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7049566063dSJacob Faibussowitsch   PetscCall(TaoAddLineSearchCounts(tao));
705c14b763aSAlp Dener 
706b2d8c577SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
707c14b763aSAlp Dener     /* Linesearch failed, revert solution */
708c14b763aSAlp Dener     bnk->f = bnk->fold;
7099566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->Xold, tao->solution));
7109566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
711c14b763aSAlp Dener 
712937a31a1SAlp Dener     switch (*stepType) {
713c14b763aSAlp Dener     case BNK_NEWTON:
7148d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
715c14b763aSAlp Dener          Update the perturbation for next time */
716c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
7172e6e4ca1SStefano Zampini         PetscBool is_gltr;
7182e6e4ca1SStefano Zampini 
719c14b763aSAlp Dener         /* Initialize the perturbation */
720c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
7219566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr));
7222e6e4ca1SStefano Zampini         if (is_gltr) {
7239566063dSJacob Faibussowitsch           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
724c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
725c14b763aSAlp Dener         }
726c14b763aSAlp Dener       } else {
727c14b763aSAlp Dener         /* Increase the perturbation */
728c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
729c14b763aSAlp Dener       }
730c14b763aSAlp Dener 
7310ad3a497SAlp Dener       if (!bnk->M) {
732c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
733c14b763aSAlp Dener            Must use gradient direction in this case */
7349566063dSJacob Faibussowitsch         PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection));
735937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
736c14b763aSAlp Dener       } else {
737c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
7389566063dSJacob Faibussowitsch         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
7398d5ead36SAlp Dener         /* Check for success (descent direction)
7408d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
7419566063dSJacob Faibussowitsch         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
7423105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
743c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
744c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
745c14b763aSAlp Dener              Use steepest descent direction (scaled) */
7469566063dSJacob Faibussowitsch           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7479566063dSJacob Faibussowitsch           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7489566063dSJacob Faibussowitsch           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
749c14b763aSAlp Dener 
750c14b763aSAlp Dener           bfgsUpdates = 1;
751937a31a1SAlp Dener           *stepType   = BNK_SCALED_GRADIENT;
752c14b763aSAlp Dener         } else {
7539566063dSJacob Faibussowitsch           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
754c14b763aSAlp Dener           if (1 == bfgsUpdates) {
755c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
756937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
757c14b763aSAlp Dener           } else {
758937a31a1SAlp Dener             *stepType = BNK_BFGS;
759c14b763aSAlp Dener           }
760c14b763aSAlp Dener         }
761c14b763aSAlp Dener       }
762c14b763aSAlp Dener       break;
763c14b763aSAlp Dener 
764c14b763aSAlp Dener     case BNK_BFGS:
765c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
766c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
767c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
7689566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7699566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7709566063dSJacob Faibussowitsch       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
771c14b763aSAlp Dener 
772c14b763aSAlp Dener       bfgsUpdates = 1;
773937a31a1SAlp Dener       *stepType   = BNK_SCALED_GRADIENT;
774c14b763aSAlp Dener       break;
775c14b763aSAlp Dener     }
7768d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7779566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
7789566063dSJacob Faibussowitsch     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
779c14b763aSAlp Dener 
7808d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
7819566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7829566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
783c14b763aSAlp Dener   }
784c14b763aSAlp Dener   *reason = ls_reason;
785c14b763aSAlp Dener   PetscFunctionReturn(0);
786c14b763aSAlp Dener }
787c14b763aSAlp Dener 
788df278d8fSAlp Dener /*------------------------------------------------------------*/
789df278d8fSAlp Dener 
790df278d8fSAlp Dener /* Routine for updating the trust radius.
791df278d8fSAlp Dener 
792df278d8fSAlp Dener   Function features three different update methods:
793df278d8fSAlp Dener   1) Line-search step length based
794df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
795df278d8fSAlp Dener   3) Interpolation
796df278d8fSAlp Dener */
797df278d8fSAlp Dener 
7989371c9d4SSatish Balay PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept) {
799080d2917SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
800080d2917SAlp Dener 
801b1c2d0e3SAlp Dener   PetscReal step, kappa;
802080d2917SAlp Dener   PetscReal gdx, tau_1, tau_2, tau_min, tau_max;
803080d2917SAlp Dener 
804080d2917SAlp Dener   PetscFunctionBegin;
805080d2917SAlp Dener   /* Update trust region radius */
806080d2917SAlp Dener   *accept = PETSC_FALSE;
80728017e9fSAlp Dener   switch (updateType) {
808080d2917SAlp Dener   case BNK_UPDATE_STEP:
809c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
810080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
8119566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step));
812080d2917SAlp Dener       if (step < bnk->nu1) {
813080d2917SAlp Dener         /* Very bad step taken; reduce radius */
814080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
815080d2917SAlp Dener       } else if (step < bnk->nu2) {
816080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
817080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
818080d2917SAlp Dener       } else if (step < bnk->nu3) {
819080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
820080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
821080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
822080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
823080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
824080d2917SAlp Dener         }
825080d2917SAlp Dener       } else if (step < bnk->nu4) {
826080d2917SAlp Dener         /*  Full step taken; increase the radius */
827080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
828080d2917SAlp Dener       } else {
829080d2917SAlp Dener         /*  More than full step taken; increase the radius */
830080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
831080d2917SAlp Dener       }
832080d2917SAlp Dener     } else {
833080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
834080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
835080d2917SAlp Dener     }
836080d2917SAlp Dener     break;
837080d2917SAlp Dener 
838080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
839080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
840e0ed867bSAlp Dener       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
841fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
842fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
843fed79b8eSAlp Dener            be rejected! */
844080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8453105154fSTodd Munson       } else {
846b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
847080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
848080d2917SAlp Dener         } else {
8493105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon)) {
850080d2917SAlp Dener             kappa = 1.0;
8513105154fSTodd Munson           } else {
852080d2917SAlp Dener             kappa = actred / prered;
853080d2917SAlp Dener           }
854fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
855080d2917SAlp Dener           if (kappa < bnk->eta1) {
856fed79b8eSAlp Dener             /* Reject the step */
857080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8583105154fSTodd Munson           } else {
859fed79b8eSAlp Dener             /* Accept the step */
860c133c014SAlp Dener             *accept = PETSC_TRUE;
861c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8628d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
863080d2917SAlp Dener               if (kappa < bnk->eta2) {
864080d2917SAlp Dener                 /* Marginal bad step */
865c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
8663105154fSTodd Munson               } else if (kappa < bnk->eta3) {
867fed79b8eSAlp Dener                 /* Reasonable step */
868fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
8693105154fSTodd Munson               } else if (kappa < bnk->eta4) {
870080d2917SAlp Dener                 /* Good step */
871c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
8723105154fSTodd Munson               } else {
873080d2917SAlp Dener                 /* Very good step */
874c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
875080d2917SAlp Dener               }
876c133c014SAlp Dener             }
877080d2917SAlp Dener           }
878080d2917SAlp Dener         }
879080d2917SAlp Dener       }
880080d2917SAlp Dener     } else {
881080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
882080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
883080d2917SAlp Dener     }
884080d2917SAlp Dener     break;
885080d2917SAlp Dener 
886080d2917SAlp Dener   default:
887080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
888b1c2d0e3SAlp Dener       if (prered < 0.0) {
889080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
890080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
891080d2917SAlp Dener         /*  be rejected! */
892080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
893080d2917SAlp Dener       } else {
894b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
895080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
896080d2917SAlp Dener         } else {
897080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
898080d2917SAlp Dener             kappa = 1.0;
899080d2917SAlp Dener           } else {
900080d2917SAlp Dener             kappa = actred / prered;
901080d2917SAlp Dener           }
902080d2917SAlp Dener 
9039566063dSJacob Faibussowitsch           PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
904080d2917SAlp Dener           tau_1   = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
905080d2917SAlp Dener           tau_2   = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
906080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
907080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
908080d2917SAlp Dener 
909080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
910080d2917SAlp Dener             /*  Great agreement */
911080d2917SAlp Dener             *accept = PETSC_TRUE;
912080d2917SAlp Dener             if (tau_max < 1.0) {
913080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
914080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
915080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
916080d2917SAlp Dener             } else {
917080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
918080d2917SAlp Dener             }
919080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
920080d2917SAlp Dener             /*  Good agreement */
921080d2917SAlp Dener             *accept = PETSC_TRUE;
922080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
923080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
924080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
925080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
926080d2917SAlp Dener             } else if (tau_max < 1.0) {
927080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
928080d2917SAlp Dener             } else {
929080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
930080d2917SAlp Dener             }
931080d2917SAlp Dener           } else {
932080d2917SAlp Dener             /*  Not good agreement */
933080d2917SAlp Dener             if (tau_min > 1.0) {
934080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
935080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
936080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
937080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
938080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
939080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
940080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
941080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
942080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
943080d2917SAlp Dener             } else {
944080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
945080d2917SAlp Dener             }
946080d2917SAlp Dener           }
947080d2917SAlp Dener         }
948080d2917SAlp Dener       }
949080d2917SAlp Dener     } else {
950080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
951080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
952080d2917SAlp Dener     }
95328017e9fSAlp Dener     break;
954080d2917SAlp Dener   }
955c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
956c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
957fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
958080d2917SAlp Dener   PetscFunctionReturn(0);
959080d2917SAlp Dener }
960080d2917SAlp Dener 
961eb910715SAlp Dener /* ---------------------------------------------------------- */
962df278d8fSAlp Dener 
9639371c9d4SSatish Balay PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType) {
96462675beeSAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
96562675beeSAlp Dener 
96662675beeSAlp Dener   PetscFunctionBegin;
96762675beeSAlp Dener   switch (stepType) {
9689371c9d4SSatish Balay   case BNK_NEWTON: ++bnk->newt; break;
9699371c9d4SSatish Balay   case BNK_BFGS: ++bnk->bfgs; break;
9709371c9d4SSatish Balay   case BNK_SCALED_GRADIENT: ++bnk->sgrad; break;
9719371c9d4SSatish Balay   case BNK_GRADIENT: ++bnk->grad; break;
9729371c9d4SSatish Balay   default: break;
97362675beeSAlp Dener   }
97462675beeSAlp Dener   PetscFunctionReturn(0);
97562675beeSAlp Dener }
97662675beeSAlp Dener 
97762675beeSAlp Dener /* ---------------------------------------------------------- */
97862675beeSAlp Dener 
9799371c9d4SSatish Balay PetscErrorCode TaoSetUp_BNK(Tao tao) {
980eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
981e031d6f5SAlp Dener   PetscInt i;
982eb910715SAlp Dener 
983eb910715SAlp Dener   PetscFunctionBegin;
984*48a46eb9SPierre Jolivet   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
985*48a46eb9SPierre Jolivet   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
986*48a46eb9SPierre Jolivet   if (!bnk->W) PetscCall(VecDuplicate(tao->solution, &bnk->W));
987*48a46eb9SPierre Jolivet   if (!bnk->Xold) PetscCall(VecDuplicate(tao->solution, &bnk->Xold));
988*48a46eb9SPierre Jolivet   if (!bnk->Gold) PetscCall(VecDuplicate(tao->solution, &bnk->Gold));
989*48a46eb9SPierre Jolivet   if (!bnk->Xwork) PetscCall(VecDuplicate(tao->solution, &bnk->Xwork));
990*48a46eb9SPierre Jolivet   if (!bnk->Gwork) PetscCall(VecDuplicate(tao->solution, &bnk->Gwork));
991*48a46eb9SPierre Jolivet   if (!bnk->unprojected_gradient) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient));
992*48a46eb9SPierre Jolivet   if (!bnk->unprojected_gradient_old) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient_old));
993*48a46eb9SPierre Jolivet   if (!bnk->Diag_min) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_min));
994*48a46eb9SPierre Jolivet   if (!bnk->Diag_max) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_max));
995e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
996c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
997c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
9989566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old)));
9999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old));
100089da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
10019566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient)));
10029566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient));
1003c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
10049566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(bnk->Gold)));
10059566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg_ctx->G_old));
1006c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
10079566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(tao->gradient)));
10089566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg->gradient));
1009c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
10109566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(tao->stepdirection)));
10119566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnk->bncg->stepdirection));
1012c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
10139566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(bnk->bncg, tao->solution));
1014c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
10159566063dSJacob Faibussowitsch     PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its));
10169566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol));
10179566063dSJacob Faibussowitsch     PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin));
10189566063dSJacob Faibussowitsch     PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP));
10199566063dSJacob Faibussowitsch     PetscCall(TaoSetObjective(bnk->bncg, tao->ops->computeobjective, tao->user_objP));
10209566063dSJacob Faibussowitsch     PetscCall(TaoSetGradient(bnk->bncg, NULL, tao->ops->computegradient, tao->user_gradP));
10219566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(bnk->bncg, NULL, tao->ops->computeobjectiveandgradient, tao->user_objgradP));
10229566063dSJacob Faibussowitsch     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg)));
1023c4b75bccSAlp Dener     for (i = 0; i < tao->numbermonitors; ++i) {
10249566063dSJacob Faibussowitsch       PetscCall(TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
10259566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(tao->monitorcontext[i])));
1026e031d6f5SAlp Dener     }
1027e031d6f5SAlp Dener   }
102883c8fe1dSLisandro Dalcin   bnk->X_inactive    = NULL;
102983c8fe1dSLisandro Dalcin   bnk->G_inactive    = NULL;
103083c8fe1dSLisandro Dalcin   bnk->inactive_work = NULL;
103183c8fe1dSLisandro Dalcin   bnk->active_work   = NULL;
103283c8fe1dSLisandro Dalcin   bnk->inactive_idx  = NULL;
103383c8fe1dSLisandro Dalcin   bnk->active_idx    = NULL;
103483c8fe1dSLisandro Dalcin   bnk->active_lower  = NULL;
103583c8fe1dSLisandro Dalcin   bnk->active_upper  = NULL;
103683c8fe1dSLisandro Dalcin   bnk->active_fixed  = NULL;
103783c8fe1dSLisandro Dalcin   bnk->M             = NULL;
103883c8fe1dSLisandro Dalcin   bnk->H_inactive    = NULL;
103983c8fe1dSLisandro Dalcin   bnk->Hpre_inactive = NULL;
1040eb910715SAlp Dener   PetscFunctionReturn(0);
1041eb910715SAlp Dener }
1042eb910715SAlp Dener 
1043eb910715SAlp Dener /*------------------------------------------------------------*/
1044df278d8fSAlp Dener 
10459371c9d4SSatish Balay PetscErrorCode TaoDestroy_BNK(Tao tao) {
1046eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
1047eb910715SAlp Dener 
1048eb910715SAlp Dener   PetscFunctionBegin;
10499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->W));
10509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Xold));
10519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Gold));
10529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Xwork));
10539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Gwork));
10549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->unprojected_gradient));
10559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->unprojected_gradient_old));
10569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Diag_min));
10579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bnk->Diag_max));
10589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_lower));
10599566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_upper));
10609566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_fixed));
10619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->active_idx));
10629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bnk->inactive_idx));
10639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->Hpre_inactive));
10649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bnk->H_inactive));
10659566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&bnk->bncg));
1066a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
10679566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
1068eb910715SAlp Dener   PetscFunctionReturn(0);
1069eb910715SAlp Dener }
1070eb910715SAlp Dener 
1071eb910715SAlp Dener /*------------------------------------------------------------*/
1072df278d8fSAlp Dener 
10739371c9d4SSatish Balay PetscErrorCode TaoSetFromOptions_BNK(Tao tao, PetscOptionItems *PetscOptionsObject) {
1074eb910715SAlp Dener   TAO_BNK *bnk = (TAO_BNK *)tao->data;
1075eb910715SAlp Dener 
1076eb910715SAlp Dener   PetscFunctionBegin;
1077d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton-Krylov method for bound constrained optimization");
10789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL));
10799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL));
10809566063dSJacob 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));
10819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval, NULL));
10829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin, NULL));
10839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax, NULL));
10849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac, NULL));
10859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin, NULL));
10869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax, NULL));
10879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac, NULL));
10889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac, NULL));
10899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac, NULL));
10909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac, NULL));
10919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1, NULL));
10929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2, NULL));
10939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3, NULL));
10949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4, NULL));
10959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1, NULL));
10969566063dSJacob 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));
10979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3, NULL));
10989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4, NULL));
10999566063dSJacob 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));
11009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1, NULL));
11019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2, NULL));
11029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3, NULL));
11039566063dSJacob 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));
11049566063dSJacob 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));
11059566063dSJacob 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));
11069566063dSJacob 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));
11079566063dSJacob 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));
11089566063dSJacob 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));
11099566063dSJacob 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));
11109566063dSJacob 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));
11119566063dSJacob 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));
11129566063dSJacob 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));
11139566063dSJacob 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));
11149566063dSJacob 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));
11159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i, NULL));
11169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1, NULL));
11179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2, NULL));
11189566063dSJacob 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));
11199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2, NULL));
11209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3, NULL));
11219566063dSJacob 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));
11229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta, NULL));
11239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius, NULL));
11249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius, NULL));
11259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
11269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
11279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
11289566063dSJacob 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));
1129d0609cedSBarry Smith   PetscOptionsHeadEnd();
11308ebe3e4eSStefano Zampini 
11319566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)(tao))->prefix));
11329566063dSJacob Faibussowitsch   PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_cg_"));
11339566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(bnk->bncg));
11348ebe3e4eSStefano Zampini 
11359566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, ((PetscObject)(tao))->prefix));
11369566063dSJacob Faibussowitsch   PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_bnk_"));
11379566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
1138eb910715SAlp Dener   PetscFunctionReturn(0);
1139eb910715SAlp Dener }
1140eb910715SAlp Dener 
1141eb910715SAlp Dener /*------------------------------------------------------------*/
1142df278d8fSAlp Dener 
11439371c9d4SSatish Balay PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) {
1144eb910715SAlp Dener   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
1145eb910715SAlp Dener   PetscInt  nrejects;
1146eb910715SAlp Dener   PetscBool isascii;
1147eb910715SAlp Dener 
1148eb910715SAlp Dener   PetscFunctionBegin;
11499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1150eb910715SAlp Dener   if (isascii) {
11519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1152b9ac7092SAlp Dener     if (bnk->M) {
11539566063dSJacob Faibussowitsch       PetscCall(MatLMVMGetRejectCount(bnk->M, &nrejects));
115463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n", nrejects));
1155eb910715SAlp Dener     }
115663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its));
115763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt));
1158*48a46eb9SPierre Jolivet     if (bnk->M) PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs));
115963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad));
116063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad));
11619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n"));
116263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  atol: %" PetscInt_FMT "\n", bnk->ksp_atol));
116363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol));
116463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol));
116563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  negc: %" PetscInt_FMT "\n", bnk->ksp_negc));
116663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol));
116763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  iter: %" PetscInt_FMT "\n", bnk->ksp_iter));
116863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  othr: %" PetscInt_FMT "\n", bnk->ksp_othr));
11699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
1170eb910715SAlp Dener   }
1171eb910715SAlp Dener   PetscFunctionReturn(0);
1172eb910715SAlp Dener }
1173eb910715SAlp Dener 
1174eb910715SAlp Dener /* ---------------------------------------------------------- */
1175df278d8fSAlp Dener 
1176eb910715SAlp Dener /*MC
1177eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
117866ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1179eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1180eb910715SAlp Dener               Hk dk = -gk
11812b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
11822b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1183eb910715SAlp Dener 
1184eb910715SAlp Dener     Options Database Keys:
11859fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
11869fa2b5dcSStefano Zampini . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
11879fa2b5dcSStefano Zampini . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
11889fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
11899fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
11909fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
11919fa2b5dcSStefano Zampini . -tao_bnk_sval - (developer) Hessian perturbation starting value
11929fa2b5dcSStefano Zampini . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
11939fa2b5dcSStefano Zampini . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
11949fa2b5dcSStefano Zampini . -tao_bnk_pmin - (developer) minimum Hessian perturbation
11959fa2b5dcSStefano Zampini . -tao_bnk_pmax - (developer) aximum Hessian perturbation
11969fa2b5dcSStefano Zampini . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
11979fa2b5dcSStefano Zampini . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
11989fa2b5dcSStefano Zampini . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
11999fa2b5dcSStefano Zampini . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
12009fa2b5dcSStefano Zampini . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
12019fa2b5dcSStefano Zampini . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction)
12029fa2b5dcSStefano Zampini . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
12039fa2b5dcSStefano Zampini . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
12049fa2b5dcSStefano Zampini . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction)
12059fa2b5dcSStefano Zampini . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
12069fa2b5dcSStefano Zampini . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
12079fa2b5dcSStefano Zampini . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
12089fa2b5dcSStefano Zampini . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
12099fa2b5dcSStefano Zampini . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
12109fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
12119fa2b5dcSStefano Zampini . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
12129fa2b5dcSStefano Zampini . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation)
12139fa2b5dcSStefano Zampini . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
12149fa2b5dcSStefano Zampini . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
12159fa2b5dcSStefano Zampini . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
12169fa2b5dcSStefano Zampini . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
12179fa2b5dcSStefano Zampini . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation)
12189fa2b5dcSStefano Zampini . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step)
12199fa2b5dcSStefano Zampini . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
12209fa2b5dcSStefano Zampini . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step)
12219fa2b5dcSStefano Zampini . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step)
12229fa2b5dcSStefano Zampini . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
12239fa2b5dcSStefano Zampini . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
12249fa2b5dcSStefano Zampini . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step)
12259fa2b5dcSStefano Zampini . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
12269fa2b5dcSStefano Zampini . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
12279fa2b5dcSStefano Zampini . -tao_bnk_mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
12289fa2b5dcSStefano Zampini . -tao_bnk_mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
12299fa2b5dcSStefano Zampini . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
12309fa2b5dcSStefano Zampini . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
12319fa2b5dcSStefano Zampini . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
12329fa2b5dcSStefano Zampini . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
12339fa2b5dcSStefano Zampini - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1234eb910715SAlp Dener 
1235eb910715SAlp Dener   Level: beginner
1236eb910715SAlp Dener M*/
1237eb910715SAlp Dener 
12389371c9d4SSatish Balay PetscErrorCode TaoCreate_BNK(Tao tao) {
1239eb910715SAlp Dener   TAO_BNK *bnk;
1240b9ac7092SAlp Dener   PC       pc;
1241eb910715SAlp Dener 
1242eb910715SAlp Dener   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao, &bnk));
1244eb910715SAlp Dener 
1245eb910715SAlp Dener   tao->ops->setup          = TaoSetUp_BNK;
1246eb910715SAlp Dener   tao->ops->view           = TaoView_BNK;
1247eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1248eb910715SAlp Dener   tao->ops->destroy        = TaoDestroy_BNK;
1249eb910715SAlp Dener 
1250eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1251eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1252eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1253eb910715SAlp Dener 
1254eb910715SAlp Dener   tao->data = (void *)bnk;
1255eb910715SAlp Dener 
125666ed3702SAlp Dener   /*  Hessian shifting parameters */
1257e0ed867bSAlp Dener   bnk->computehessian = TaoBNKComputeHessian;
1258e0ed867bSAlp Dener   bnk->computestep    = TaoBNKComputeStep;
1259e0ed867bSAlp Dener 
1260eb910715SAlp Dener   bnk->sval  = 0.0;
1261eb910715SAlp Dener   bnk->imin  = 1.0e-4;
1262eb910715SAlp Dener   bnk->imax  = 1.0e+2;
1263eb910715SAlp Dener   bnk->imfac = 1.0e-1;
1264eb910715SAlp Dener 
1265eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1266eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1267eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1268eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1269eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1270eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1271eb910715SAlp Dener 
1272eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1273eb910715SAlp Dener   bnk->nu1 = 0.25;
1274eb910715SAlp Dener   bnk->nu2 = 0.50;
1275eb910715SAlp Dener   bnk->nu3 = 1.00;
1276eb910715SAlp Dener   bnk->nu4 = 1.25;
1277eb910715SAlp Dener 
1278eb910715SAlp Dener   bnk->omega1 = 0.25;
1279eb910715SAlp Dener   bnk->omega2 = 0.50;
1280eb910715SAlp Dener   bnk->omega3 = 1.00;
1281eb910715SAlp Dener   bnk->omega4 = 2.00;
1282eb910715SAlp Dener   bnk->omega5 = 4.00;
1283eb910715SAlp Dener 
1284eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1285eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1286eb910715SAlp Dener   bnk->eta2 = 0.25;
1287eb910715SAlp Dener   bnk->eta3 = 0.50;
1288eb910715SAlp Dener   bnk->eta4 = 0.90;
1289eb910715SAlp Dener 
1290eb910715SAlp Dener   bnk->alpha1 = 0.25;
1291eb910715SAlp Dener   bnk->alpha2 = 0.50;
1292eb910715SAlp Dener   bnk->alpha3 = 1.00;
1293eb910715SAlp Dener   bnk->alpha4 = 2.00;
1294eb910715SAlp Dener   bnk->alpha5 = 4.00;
1295eb910715SAlp Dener 
1296eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1297eb910715SAlp Dener   bnk->mu1 = 0.10;
1298eb910715SAlp Dener   bnk->mu2 = 0.50;
1299eb910715SAlp Dener 
1300eb910715SAlp Dener   bnk->gamma1 = 0.25;
1301eb910715SAlp Dener   bnk->gamma2 = 0.50;
1302eb910715SAlp Dener   bnk->gamma3 = 2.00;
1303eb910715SAlp Dener   bnk->gamma4 = 4.00;
1304eb910715SAlp Dener 
1305eb910715SAlp Dener   bnk->theta = 0.05;
1306eb910715SAlp Dener 
1307eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1308eb910715SAlp Dener   bnk->mu1_i = 0.35;
1309eb910715SAlp Dener   bnk->mu2_i = 0.50;
1310eb910715SAlp Dener 
1311eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1312eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1313eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1314eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1315eb910715SAlp Dener 
1316eb910715SAlp Dener   bnk->theta_i = 0.25;
1317eb910715SAlp Dener 
1318eb910715SAlp Dener   /*  Remaining parameters */
1319c0f10754SAlp Dener   bnk->max_cg_its = 0;
1320eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1321eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1322770b7498SAlp Dener   bnk->epsilon    = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0);
13230a4511e9SAlp Dener   bnk->as_tol     = 1.0e-3;
13240a4511e9SAlp Dener   bnk->as_step    = 1.0e-3;
132562675beeSAlp Dener   bnk->dmin       = 1.0e-6;
132662675beeSAlp Dener   bnk->dmax       = 1.0e6;
1327eb910715SAlp Dener 
132883c8fe1dSLisandro Dalcin   bnk->M           = NULL;
132983c8fe1dSLisandro Dalcin   bnk->bfgs_pre    = NULL;
1330eb910715SAlp Dener   bnk->init_type   = BNK_INIT_INTERPOLATION;
13317b1c7716SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION;
13322f75a4aaSAlp Dener   bnk->as_type     = BNK_AS_BERTSEKAS;
1333eb910715SAlp Dener 
1334e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
13359566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg));
13369566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1));
13379566063dSJacob Faibussowitsch   PetscCall(TaoSetType(bnk->bncg, TAOBNCG));
1338e031d6f5SAlp Dener 
1339c0f10754SAlp Dener   /* Create the line search */
13409566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
13419566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1342f4db9bf7SStefano Zampini   PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
13439566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
1344eb910715SAlp Dener 
1345eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
13469566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
13479566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
13489566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
13499566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(tao->ksp, &pc));
13509566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLMVM));
1351eb910715SAlp Dener   PetscFunctionReturn(0);
1352eb910715SAlp Dener }
1353