1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener #include <petscksp.h>
4eb910715SAlp Dener
570a3f44bSAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
670a3f44bSAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
770a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
870a3f44bSAlp Dener
9b3e6a353SBarry Smith /* Extracts from the full Hessian the part associated with the current bnk->inactive_idx and set the PCLMVM preconditioner */
10e031d6f5SAlp Dener
TaoBNKComputeSubHessian(Tao tao)11b3e6a353SBarry Smith static PetscErrorCode TaoBNKComputeSubHessian(Tao tao)
12b3e6a353SBarry Smith {
13b3e6a353SBarry Smith TAO_BNK *bnk = (TAO_BNK *)tao->data;
14b3e6a353SBarry Smith
15b3e6a353SBarry Smith PetscFunctionBegin;
16b3e6a353SBarry Smith PetscCall(MatDestroy(&bnk->Hpre_inactive));
17b3e6a353SBarry Smith PetscCall(MatDestroy(&bnk->H_inactive));
18b3e6a353SBarry Smith if (bnk->active_idx) {
19b3e6a353SBarry Smith PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
20b3e6a353SBarry Smith if (tao->hessian == tao->hessian_pre) {
21b3e6a353SBarry Smith PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
22b3e6a353SBarry Smith bnk->Hpre_inactive = bnk->H_inactive;
23b3e6a353SBarry Smith } else {
24b3e6a353SBarry Smith PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive));
25b3e6a353SBarry Smith }
26b3e6a353SBarry Smith if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx));
27b3e6a353SBarry Smith } else {
28b3e6a353SBarry Smith PetscCall(PetscObjectReference((PetscObject)tao->hessian));
29b3e6a353SBarry Smith bnk->H_inactive = tao->hessian;
30b3e6a353SBarry Smith PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre));
31b3e6a353SBarry Smith bnk->Hpre_inactive = tao->hessian_pre;
32b3e6a353SBarry Smith if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre));
33b3e6a353SBarry Smith }
343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
35b3e6a353SBarry Smith }
36b3e6a353SBarry Smith
37b3e6a353SBarry Smith /* Initializes the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
38df278d8fSAlp Dener
TaoBNKInitialize(Tao tao,PetscInt initType,PetscBool * needH)39d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
40d71ae5a4SJacob Faibussowitsch {
41eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
42eb910715SAlp Dener PC pc;
4389da521bSAlp Dener PetscReal f_min, ftrial, prered, actred, kappa, sigma, resnorm;
44eb910715SAlp Dener PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
450ad3a497SAlp Dener PetscBool is_bfgs, is_jacobi, is_symmetric, sym_set;
46c4b75bccSAlp Dener PetscInt n, N, nDiff;
47eb910715SAlp Dener PetscInt i_max = 5;
48eb910715SAlp Dener PetscInt j_max = 1;
49eb910715SAlp Dener PetscInt i, j;
500cd8b6e2SStefano Zampini PetscBool kspTR;
51eb910715SAlp Dener
52eb910715SAlp Dener PetscFunctionBegin;
5328017e9fSAlp Dener /* Project the current point onto the feasible set */
549566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao));
559566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU));
561baa6e33SBarry Smith if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
5728017e9fSAlp Dener
5828017e9fSAlp Dener /* Project the initial point onto the feasible region */
599566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
6028017e9fSAlp Dener
6128017e9fSAlp Dener /* Check convergence criteria */
629566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
639566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
649566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
65976ed0a4SStefano Zampini if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
669566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
6728017e9fSAlp Dener
68c0f10754SAlp Dener /* Test the initial point for convergence */
699566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
709566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
71*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
729566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
739566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
74dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
753ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
76c0f10754SAlp Dener
77e031d6f5SAlp Dener /* Reset KSP stopping reason counters */
78eb910715SAlp Dener bnk->ksp_atol = 0;
79eb910715SAlp Dener bnk->ksp_rtol = 0;
80eb910715SAlp Dener bnk->ksp_dtol = 0;
81eb910715SAlp Dener bnk->ksp_ctol = 0;
82eb910715SAlp Dener bnk->ksp_negc = 0;
83eb910715SAlp Dener bnk->ksp_iter = 0;
84eb910715SAlp Dener bnk->ksp_othr = 0;
85eb910715SAlp Dener
86e031d6f5SAlp Dener /* Reset accepted step type counters */
87e031d6f5SAlp Dener bnk->tot_cg_its = 0;
88e031d6f5SAlp Dener bnk->newt = 0;
89e031d6f5SAlp Dener bnk->bfgs = 0;
90e031d6f5SAlp Dener bnk->sgrad = 0;
91e031d6f5SAlp Dener bnk->grad = 0;
92e031d6f5SAlp Dener
93fed79b8eSAlp Dener /* Initialize the Hessian perturbation */
94fed79b8eSAlp Dener bnk->pert = bnk->sval;
95fed79b8eSAlp Dener
96937a31a1SAlp Dener /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
979566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0));
98937a31a1SAlp Dener
99e031d6f5SAlp Dener /* Allocate the vectors needed for the BFGS approximation */
1009566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc));
1019566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
1029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
103b9ac7092SAlp Dener if (is_bfgs) {
104b9ac7092SAlp Dener bnk->bfgs_pre = pc;
1059566063dSJacob Faibussowitsch PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M));
1069566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->solution, &n));
1079566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &N));
1089566063dSJacob Faibussowitsch PetscCall(MatSetSizes(bnk->M, n, n, N, N));
1099566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient));
1109566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric));
1113c859ba3SBarry Smith PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
1121baa6e33SBarry Smith } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));
113e031d6f5SAlp Dener
114e031d6f5SAlp Dener /* Prepare the min/max vectors for safeguarding diagonal scales */
1159566063dSJacob Faibussowitsch PetscCall(VecSet(bnk->Diag_min, bnk->dmin));
1169566063dSJacob Faibussowitsch PetscCall(VecSet(bnk->Diag_max, bnk->dmax));
117eb910715SAlp Dener
118eb910715SAlp Dener /* Initialize trust-region radius. The initialization is only performed
119eb910715SAlp Dener when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
120c0f10754SAlp Dener *needH = PETSC_TRUE;
1210cd8b6e2SStefano Zampini PetscCall(PetscObjectHasFunction((PetscObject)tao->ksp, "KSPCGSetRadius_C", &kspTR));
1222e6e4ca1SStefano Zampini if (kspTR) {
12362675beeSAlp Dener switch (initType) {
124eb910715SAlp Dener case BNK_INIT_CONSTANT:
125eb910715SAlp Dener /* Use the initial radius specified */
126c0f10754SAlp Dener tao->trust = tao->trust0;
127eb910715SAlp Dener break;
128eb910715SAlp Dener
129eb910715SAlp Dener case BNK_INIT_INTERPOLATION:
130c0f10754SAlp Dener /* Use interpolation based on the initial Hessian */
131eb910715SAlp Dener max_radius = 0.0;
13208752603SAlp Dener tao->trust = tao->trust0;
133eb910715SAlp Dener for (j = 0; j < j_max; ++j) {
1340a4511e9SAlp Dener f_min = bnk->f;
135eb910715SAlp Dener sigma = 0.0;
136eb910715SAlp Dener
137c0f10754SAlp Dener if (*needH) {
13862602cfbSAlp Dener /* Compute the Hessian at the new step, and extract the inactive subsystem */
1399566063dSJacob Faibussowitsch PetscCall((*bnk->computehessian)(tao));
1409566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE));
141b3e6a353SBarry Smith PetscCall(TaoBNKComputeSubHessian(tao));
142c0f10754SAlp Dener *needH = PETSC_FALSE;
143eb910715SAlp Dener }
144eb910715SAlp Dener
145eb910715SAlp Dener for (i = 0; i < i_max; ++i) {
14662602cfbSAlp Dener /* Take a steepest descent step and snap it to bounds */
1479566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold));
1489566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, -tao->trust / bnk->gnorm, tao->gradient));
1499566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
15089da521bSAlp Dener /* Compute the step we actually accepted */
1519566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->W));
1529566063dSJacob Faibussowitsch PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold));
15362602cfbSAlp Dener /* Compute the objective at the trial */
1549566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial));
155*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
1569566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution));
157eb910715SAlp Dener if (PetscIsInfOrNanReal(ftrial)) {
158eb910715SAlp Dener tau = bnk->gamma1_i;
159eb910715SAlp Dener } else {
1600a4511e9SAlp Dener if (ftrial < f_min) {
1610a4511e9SAlp Dener f_min = ftrial;
162eb910715SAlp Dener sigma = -tao->trust / bnk->gnorm;
163eb910715SAlp Dener }
16408752603SAlp Dener
165770b7498SAlp Dener /* Compute the predicted and actual reduction */
16689da521bSAlp Dener if (bnk->active_idx) {
1679566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1689566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1692ab2a32cSAlp Dener } else {
17008752603SAlp Dener bnk->X_inactive = bnk->W;
17108752603SAlp Dener bnk->inactive_work = bnk->Xwork;
1722ab2a32cSAlp Dener }
1739566063dSJacob Faibussowitsch PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
1749566063dSJacob Faibussowitsch PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered));
17589da521bSAlp Dener if (bnk->active_idx) {
1769566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
1779566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
1782ab2a32cSAlp Dener }
179eb910715SAlp Dener prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
180eb910715SAlp Dener actred = bnk->f - ftrial;
1813105154fSTodd Munson if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
182eb910715SAlp Dener kappa = 1.0;
1833105154fSTodd Munson } else {
184eb910715SAlp Dener kappa = actred / prered;
185eb910715SAlp Dener }
186eb910715SAlp Dener
187eb910715SAlp Dener tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
188eb910715SAlp Dener tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
189eb910715SAlp Dener tau_min = PetscMin(tau_1, tau_2);
190eb910715SAlp Dener tau_max = PetscMax(tau_1, tau_2);
191eb910715SAlp Dener
19218cfbf8eSSatish Balay if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
193eb910715SAlp Dener /* Great agreement */
194eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust);
195eb910715SAlp Dener
196eb910715SAlp Dener if (tau_max < 1.0) {
197eb910715SAlp Dener tau = bnk->gamma3_i;
1983105154fSTodd Munson } else if (tau_max > bnk->gamma4_i) {
199eb910715SAlp Dener tau = bnk->gamma4_i;
2003105154fSTodd Munson } else {
201eb910715SAlp Dener tau = tau_max;
202eb910715SAlp Dener }
20318cfbf8eSSatish Balay } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
204eb910715SAlp Dener /* Good agreement */
205eb910715SAlp Dener max_radius = PetscMax(max_radius, tao->trust);
206eb910715SAlp Dener
207eb910715SAlp Dener if (tau_max < bnk->gamma2_i) {
208eb910715SAlp Dener tau = bnk->gamma2_i;
209eb910715SAlp Dener } else if (tau_max > bnk->gamma3_i) {
210eb910715SAlp Dener tau = bnk->gamma3_i;
211eb910715SAlp Dener } else {
212eb910715SAlp Dener tau = tau_max;
213eb910715SAlp Dener }
2148f8a4e06SAlp Dener } else {
215eb910715SAlp Dener /* Not good agreement */
216eb910715SAlp Dener if (tau_min > 1.0) {
217eb910715SAlp Dener tau = bnk->gamma2_i;
218eb910715SAlp Dener } else if (tau_max < bnk->gamma1_i) {
219eb910715SAlp Dener tau = bnk->gamma1_i;
220eb910715SAlp Dener } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
221eb910715SAlp Dener tau = bnk->gamma1_i;
2223105154fSTodd Munson } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
223eb910715SAlp Dener tau = tau_1;
2243105154fSTodd Munson } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
225eb910715SAlp Dener tau = tau_2;
226eb910715SAlp Dener } else {
227eb910715SAlp Dener tau = tau_max;
228eb910715SAlp Dener }
229eb910715SAlp Dener }
230eb910715SAlp Dener }
231eb910715SAlp Dener tao->trust = tau * tao->trust;
232eb910715SAlp Dener }
233eb910715SAlp Dener
2340a4511e9SAlp Dener if (f_min < bnk->f) {
235937a31a1SAlp Dener /* We accidentally found a solution better than the initial, so accept it */
2360a4511e9SAlp Dener bnk->f = f_min;
2379566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, bnk->Xold));
2389566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
2399566063dSJacob Faibussowitsch PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
2409566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tao->stepdirection));
2419566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
2429566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
2439566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
2449566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
245976ed0a4SStefano Zampini if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
246937a31a1SAlp Dener /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
2479566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
248c0f10754SAlp Dener *needH = PETSC_TRUE;
249937a31a1SAlp Dener /* Test the new step for convergence */
2509566063dSJacob Faibussowitsch PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
2519566063dSJacob Faibussowitsch PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
252*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
2539566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
2549566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
255dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
2563ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
257937a31a1SAlp Dener /* active BNCG recycling early because we have a stepdirection computed */
2589566063dSJacob Faibussowitsch PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
259eb910715SAlp Dener }
260eb910715SAlp Dener }
261eb910715SAlp Dener tao->trust = PetscMax(tao->trust, max_radius);
262e031d6f5SAlp Dener
263e031d6f5SAlp Dener /* Ensure that the trust radius is within the limits */
264e031d6f5SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius);
265e031d6f5SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius);
266eb910715SAlp Dener break;
267eb910715SAlp Dener
268eb910715SAlp Dener default:
269eb910715SAlp Dener /* Norm of the first direction will initialize radius */
270eb910715SAlp Dener tao->trust = 0.0;
271eb910715SAlp Dener break;
272eb910715SAlp Dener }
273eb910715SAlp Dener }
2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
275eb910715SAlp Dener }
276eb910715SAlp Dener
277b3e6a353SBarry Smith /* Computes the exact Hessian and extracts its subHessian */
27862675beeSAlp Dener
TaoBNKComputeHessian(Tao tao)279d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKComputeHessian(Tao tao)
280d71ae5a4SJacob Faibussowitsch {
28162675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
28262675beeSAlp Dener
28362675beeSAlp Dener PetscFunctionBegin;
28462675beeSAlp Dener /* Compute the Hessian */
2859566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
28662675beeSAlp Dener /* Add a correction to the BFGS preconditioner */
2871baa6e33SBarry Smith if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
288e0ed867bSAlp Dener /* Prepare the reduced sub-matrices for the inactive set */
289b3e6a353SBarry Smith PetscCall(TaoBNKComputeSubHessian(tao));
2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29162675beeSAlp Dener }
29262675beeSAlp Dener
2932f75a4aaSAlp Dener /* Routine for estimating the active set */
2942f75a4aaSAlp Dener
TaoBNKEstimateActiveSet(Tao tao,PetscInt asType)295d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
296d71ae5a4SJacob Faibussowitsch {
2972f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
298f4db9bf7SStefano Zampini PetscBool hessComputed, diagExists, hadactive;
2992f75a4aaSAlp Dener
3002f75a4aaSAlp Dener PetscFunctionBegin;
301f4db9bf7SStefano Zampini hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE;
30208752603SAlp Dener switch (asType) {
3032f75a4aaSAlp Dener case BNK_AS_NONE:
3049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->inactive_idx));
3059566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx));
3069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_idx));
3079566063dSJacob Faibussowitsch PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx));
3082f75a4aaSAlp Dener break;
3092f75a4aaSAlp Dener
3102f75a4aaSAlp Dener case BNK_AS_BERTSEKAS:
3112f75a4aaSAlp Dener /* Compute the trial step vector with which we will estimate the active set at the next iteration */
312b9ac7092SAlp Dener if (bnk->M) {
3137addb90fSBarry Smith /* If the BFGS matrix is available, we will construct a trial step with it */
3149566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W));
3152f75a4aaSAlp Dener } else {
316fc5ca067SStefano Zampini hessComputed = diagExists = PETSC_FALSE;
31748a46eb9SPierre Jolivet if (tao->hessian) PetscCall(MatAssembled(tao->hessian, &hessComputed));
31848a46eb9SPierre Jolivet if (hessComputed) PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists));
319fc5ca067SStefano Zampini if (diagExists) {
3209b6ef848SAlp Dener /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3219566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork));
3229566063dSJacob Faibussowitsch PetscCall(VecAbs(bnk->Xwork));
3239566063dSJacob Faibussowitsch PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork));
3249566063dSJacob Faibussowitsch PetscCall(VecReciprocal(bnk->Xwork));
3259566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient));
32661be54a6SAlp Dener } else {
327c4b75bccSAlp Dener /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
3289566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W));
32961be54a6SAlp Dener }
3302f75a4aaSAlp Dener }
3319566063dSJacob Faibussowitsch PetscCall(VecScale(bnk->W, -1.0));
3329371c9d4SSatish 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));
333c4b75bccSAlp Dener break;
3342f75a4aaSAlp Dener
335d71ae5a4SJacob Faibussowitsch default:
336d71ae5a4SJacob Faibussowitsch break;
3372f75a4aaSAlp Dener }
338f4db9bf7SStefano Zampini bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */
3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3402f75a4aaSAlp Dener }
3412f75a4aaSAlp Dener
3422f75a4aaSAlp Dener /* Routine for bounding the step direction */
3432f75a4aaSAlp Dener
TaoBNKBoundStep(Tao tao,PetscInt asType,Vec step)344d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
345d71ae5a4SJacob Faibussowitsch {
3462f75a4aaSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
3472f75a4aaSAlp Dener
3482f75a4aaSAlp Dener PetscFunctionBegin;
349a1318120SAlp Dener switch (asType) {
350d71ae5a4SJacob Faibussowitsch case BNK_AS_NONE:
351976ed0a4SStefano Zampini if (bnk->active_idx) PetscCall(VecISSet(step, bnk->active_idx, 0.0));
352d71ae5a4SJacob Faibussowitsch break;
353d71ae5a4SJacob Faibussowitsch case BNK_AS_BERTSEKAS:
354d71ae5a4SJacob Faibussowitsch PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step));
355d71ae5a4SJacob Faibussowitsch break;
356d71ae5a4SJacob Faibussowitsch default:
357d71ae5a4SJacob Faibussowitsch break;
3582f75a4aaSAlp Dener }
3593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3602f75a4aaSAlp Dener }
3612f75a4aaSAlp Dener
362e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
363e031d6f5SAlp Dener accelerate Newton convergence.
364e031d6f5SAlp Dener
365e031d6f5SAlp Dener In practice, this approach simply trades off Hessian evaluations
366e031d6f5SAlp Dener for more gradient evaluations.
367e031d6f5SAlp Dener */
368e031d6f5SAlp Dener
TaoBNKTakeCGSteps(Tao tao,PetscBool * terminate)369d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
370d71ae5a4SJacob Faibussowitsch {
371c0f10754SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
372c0f10754SAlp Dener
373c0f10754SAlp Dener PetscFunctionBegin;
374c0f10754SAlp Dener *terminate = PETSC_FALSE;
375c0f10754SAlp Dener if (bnk->max_cg_its > 0) {
376c4b75bccSAlp Dener /* Copy the current function value (important vectors are already shared) */
377c0f10754SAlp Dener bnk->bncg_ctx->f = bnk->f;
378c0f10754SAlp Dener /* Take some small finite number of BNCG iterations */
3799566063dSJacob Faibussowitsch PetscCall(TaoSolve(bnk->bncg));
380c0f10754SAlp Dener /* Add the number of gradient and function evaluations to the total */
381c0f10754SAlp Dener tao->nfuncs += bnk->bncg->nfuncs;
382c0f10754SAlp Dener tao->nfuncgrads += bnk->bncg->nfuncgrads;
383c0f10754SAlp Dener tao->ngrads += bnk->bncg->ngrads;
384c0f10754SAlp Dener tao->nhess += bnk->bncg->nhess;
385e031d6f5SAlp Dener bnk->tot_cg_its += bnk->bncg->niter;
386c4b75bccSAlp Dener /* Extract the BNCG function value out and save it into BNK */
387c0f10754SAlp Dener bnk->f = bnk->bncg_ctx->f;
388c0f10754SAlp 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) {
389c0f10754SAlp Dener *terminate = PETSC_TRUE;
39061be54a6SAlp Dener } else {
3919566063dSJacob Faibussowitsch PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
392c0f10754SAlp Dener }
393c0f10754SAlp Dener }
3943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
395c0f10754SAlp Dener }
396c0f10754SAlp Dener
397c0f10754SAlp Dener /* Routine for computing the Newton step. */
398df278d8fSAlp Dener
TaoBNKComputeStep(Tao tao,PetscBool shift,KSPConvergedReason * ksp_reason,PetscInt * step_type)399d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
400d71ae5a4SJacob Faibussowitsch {
401eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
402eb910715SAlp Dener PetscInt bfgsUpdates = 0;
403eb910715SAlp Dener PetscInt kspits;
404bddd1ffdSAlp Dener PetscBool is_lmvm;
4050cd8b6e2SStefano Zampini PetscBool kspTR;
406eb910715SAlp Dener
407eb910715SAlp Dener PetscFunctionBegin;
40889da521bSAlp Dener /* If there are no inactive variables left, save some computation and return an adjusted zero step
40989da521bSAlp Dener that has (l-x) and (u-x) for lower and upper bounded variables. */
41089da521bSAlp Dener if (!bnk->inactive_idx) {
4119566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0));
4129566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
4133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
41489da521bSAlp Dener }
41589da521bSAlp Dener
41662675beeSAlp Dener /* Shift the reduced Hessian matrix */
417e831869dSStefano Zampini if (shift && bnk->pert > 0) {
4189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm));
419f7bf01afSAlp Dener if (is_lmvm) {
4209566063dSJacob Faibussowitsch PetscCall(MatShift(tao->hessian, bnk->pert));
421f7bf01afSAlp Dener } else {
4229566063dSJacob Faibussowitsch PetscCall(MatShift(bnk->H_inactive, bnk->pert));
42348a46eb9SPierre Jolivet if (bnk->H_inactive != bnk->Hpre_inactive) PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert));
42462675beeSAlp Dener }
425f7bf01afSAlp Dener }
42662675beeSAlp Dener
427eb910715SAlp Dener /* Solve the Newton system of equations */
428937a31a1SAlp Dener tao->ksp_its = 0;
4299566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0));
430f4db9bf7SStefano Zampini if (bnk->resetksp) {
4319566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp));
4329566063dSJacob Faibussowitsch PetscCall(KSPResetFromOptions(tao->ksp));
433f4db9bf7SStefano Zampini bnk->resetksp = PETSC_FALSE;
434f4db9bf7SStefano Zampini }
4359566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, bnk->H_inactive, bnk->Hpre_inactive));
4369566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork));
43789da521bSAlp Dener if (bnk->active_idx) {
4389566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4399566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4405e9b73cbSAlp Dener } else {
4415e9b73cbSAlp Dener bnk->G_inactive = bnk->unprojected_gradient;
4425e9b73cbSAlp Dener bnk->X_inactive = tao->stepdirection;
44328017e9fSAlp Dener }
4449566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
4459566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4469566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
447eb910715SAlp Dener tao->ksp_its += kspits;
448eb910715SAlp Dener tao->ksp_tot_its += kspits;
4490cd8b6e2SStefano Zampini PetscCall(PetscObjectHasFunction((PetscObject)tao->ksp, "KSPCGGetNormD_C", &kspTR));
450f4db9bf7SStefano Zampini if (kspTR) {
4519566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
452eb910715SAlp Dener
453eb910715SAlp Dener if (0.0 == tao->trust) {
454eb910715SAlp Dener /* Radius was uninitialized; use the norm of the direction */
455080d2917SAlp Dener if (bnk->dnorm > 0.0) {
456080d2917SAlp Dener tao->trust = bnk->dnorm;
457eb910715SAlp Dener
458eb910715SAlp Dener /* Modify the radius if it is too large or small */
459eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius);
460eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius);
461eb910715SAlp Dener } else {
462eb910715SAlp Dener /* The direction was bad; set radius to default value and re-solve
463eb910715SAlp Dener the trust-region subproblem to get a direction */
464eb910715SAlp Dener tao->trust = tao->trust0;
465eb910715SAlp Dener
466eb910715SAlp Dener /* Modify the radius if it is too large or small */
467eb910715SAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius);
468eb910715SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius);
469eb910715SAlp Dener
4709566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
4719566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
4729566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
473eb910715SAlp Dener tao->ksp_its += kspits;
474eb910715SAlp Dener tao->ksp_tot_its += kspits;
4759566063dSJacob Faibussowitsch PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));
476eb910715SAlp Dener
4773c859ba3SBarry Smith PetscCheck(bnk->dnorm != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
478eb910715SAlp Dener }
479eb910715SAlp Dener }
480eb910715SAlp Dener }
4815e9b73cbSAlp Dener /* Restore sub vectors back */
48289da521bSAlp Dener if (bnk->active_idx) {
4839566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
4849566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
4855e9b73cbSAlp Dener }
486770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
4879566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
4889566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
489770b7498SAlp Dener
490770b7498SAlp Dener /* Record convergence reasons */
4919566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason));
492e465cd6fSAlp Dener if (KSP_CONVERGED_ATOL == *ksp_reason) {
493770b7498SAlp Dener ++bnk->ksp_atol;
494e465cd6fSAlp Dener } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
495770b7498SAlp Dener ++bnk->ksp_rtol;
4964a221d59SStefano Zampini } else if (KSP_CONVERGED_STEP_LENGTH == *ksp_reason) {
497770b7498SAlp Dener ++bnk->ksp_ctol;
4984a221d59SStefano Zampini } else if (KSP_CONVERGED_NEG_CURVE == *ksp_reason) {
499770b7498SAlp Dener ++bnk->ksp_negc;
500e465cd6fSAlp Dener } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
501770b7498SAlp Dener ++bnk->ksp_dtol;
502e465cd6fSAlp Dener } else if (KSP_DIVERGED_ITS == *ksp_reason) {
503770b7498SAlp Dener ++bnk->ksp_iter;
504770b7498SAlp Dener } else {
505770b7498SAlp Dener ++bnk->ksp_othr;
506770b7498SAlp Dener }
507fed79b8eSAlp Dener
508fed79b8eSAlp Dener /* Make sure the BFGS preconditioner is healthy */
509b9ac7092SAlp Dener if (bnk->M) {
5109566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
511b2d8c577SAlp Dener if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
512fed79b8eSAlp Dener /* Preconditioner is numerically indefinite; reset the approximation. */
5139566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
5149566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
515eb910715SAlp Dener }
516fed79b8eSAlp Dener }
5176b591159SAlp Dener *step_type = BNK_NEWTON;
5183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
519e465cd6fSAlp Dener }
520eb910715SAlp Dener
5215e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
5225e9b73cbSAlp Dener
TaoBNKRecomputePred(Tao tao,Vec S,PetscReal * prered)523d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
524d71ae5a4SJacob Faibussowitsch {
5255e9b73cbSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
5265e9b73cbSAlp Dener
5275e9b73cbSAlp Dener PetscFunctionBegin;
5285e9b73cbSAlp Dener /* Extract subvectors associated with the inactive set */
52989da521bSAlp Dener if (bnk->active_idx) {
5309566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5319566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5329566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5335e9b73cbSAlp Dener } else {
5345e9b73cbSAlp Dener bnk->X_inactive = tao->stepdirection;
5355e9b73cbSAlp Dener bnk->inactive_work = bnk->Xwork;
5365e9b73cbSAlp Dener bnk->G_inactive = bnk->Gwork;
5375e9b73cbSAlp Dener }
5385e9b73cbSAlp Dener /* Recompute the predicted decrease based on the quadratic model */
5399566063dSJacob Faibussowitsch PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
5409566063dSJacob Faibussowitsch PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive));
5419566063dSJacob Faibussowitsch PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered));
5425e9b73cbSAlp Dener /* Restore the sub vectors */
54389da521bSAlp Dener if (bnk->active_idx) {
5449566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
5459566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
5469566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
5475e9b73cbSAlp Dener }
5483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5495e9b73cbSAlp Dener }
5505e9b73cbSAlp Dener
55162675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
55262675beeSAlp Dener
55362675beeSAlp Dener The step direction falls back onto BFGS, scaled gradient and gradient steps
55462675beeSAlp Dener in the event that the Newton step fails the test.
55562675beeSAlp Dener */
55662675beeSAlp Dener
TaoBNKSafeguardStep(Tao tao,KSPConvergedReason ksp_reason,PetscInt * stepType)557d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
558d71ae5a4SJacob Faibussowitsch {
559e465cd6fSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
560b2d8c577SAlp Dener PetscReal gdx, e_min;
561e465cd6fSAlp Dener PetscInt bfgsUpdates;
562e465cd6fSAlp Dener
563e465cd6fSAlp Dener PetscFunctionBegin;
5646b591159SAlp Dener switch (*stepType) {
5656b591159SAlp Dener case BNK_NEWTON:
5669566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
567eb910715SAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
568*76c63389SBarry Smith /* Newton step is not descent or direction produced infinity or NaN
569eb910715SAlp Dener Update the perturbation for next time */
570eb910715SAlp Dener if (bnk->pert <= 0.0) {
5712e6e4ca1SStefano Zampini PetscBool is_gltr;
5722e6e4ca1SStefano Zampini
573eb910715SAlp Dener /* Initialize the perturbation */
574eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
575f4f49eeaSPierre Jolivet PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
5762e6e4ca1SStefano Zampini if (is_gltr) {
5779566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
578eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min);
579eb910715SAlp Dener }
580eb910715SAlp Dener } else {
581eb910715SAlp Dener /* Increase the perturbation */
582eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
583eb910715SAlp Dener }
584eb910715SAlp Dener
5850ad3a497SAlp Dener if (!bnk->M) {
586eb910715SAlp Dener /* We don't have the bfgs matrix around and updated
587eb910715SAlp Dener Must use gradient direction in this case */
5889566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tao->stepdirection));
589eb910715SAlp Dener *stepType = BNK_GRADIENT;
590eb910715SAlp Dener } else {
591eb910715SAlp Dener /* Attempt to use the BFGS direction */
5929566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
593eb910715SAlp Dener
5948d5ead36SAlp Dener /* Check for success (descent direction)
5958d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because
5968d5ead36SAlp Dener the fall-back step is missing a negative sign. */
5979566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
5983105154fSTodd Munson if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
599eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number
600eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case because
601eb910715SAlp Dener the first solve produces the scaled gradient direction,
602eb910715SAlp Dener which is guaranteed to be descent */
603eb910715SAlp Dener
604eb910715SAlp Dener /* Use steepest descent direction (scaled) */
6059566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6069566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6079566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
608eb910715SAlp Dener
609eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT;
610eb910715SAlp Dener } else {
6119566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
612eb910715SAlp Dener if (1 == bfgsUpdates) {
613eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */
614eb910715SAlp Dener *stepType = BNK_SCALED_GRADIENT;
615eb910715SAlp Dener } else {
616eb910715SAlp Dener *stepType = BNK_BFGS;
617eb910715SAlp Dener }
618eb910715SAlp Dener }
619eb910715SAlp Dener }
6208d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6219566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
6229566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
623eb910715SAlp Dener } else {
624eb910715SAlp Dener /* Computed Newton step is descent */
625eb910715SAlp Dener switch (ksp_reason) {
626eb910715SAlp Dener case KSP_DIVERGED_NANORINF:
627eb910715SAlp Dener case KSP_DIVERGED_BREAKDOWN:
628eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_MAT:
629eb910715SAlp Dener case KSP_DIVERGED_INDEFINITE_PC:
6304a221d59SStefano Zampini case KSP_CONVERGED_NEG_CURVE:
631eb910715SAlp Dener /* Matrix or preconditioner is indefinite; increase perturbation */
632eb910715SAlp Dener if (bnk->pert <= 0.0) {
6332e6e4ca1SStefano Zampini PetscBool is_gltr;
6342e6e4ca1SStefano Zampini
635eb910715SAlp Dener /* Initialize the perturbation */
636eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
637f4f49eeaSPierre Jolivet PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
6382e6e4ca1SStefano Zampini if (is_gltr) {
6399566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
640eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min);
641eb910715SAlp Dener }
642eb910715SAlp Dener } else {
643eb910715SAlp Dener /* Increase the perturbation */
644eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
645eb910715SAlp Dener }
646eb910715SAlp Dener break;
647eb910715SAlp Dener
648eb910715SAlp Dener default:
649eb910715SAlp Dener /* Newton step computation is good; decrease perturbation */
650eb910715SAlp Dener bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
651ad540459SPierre Jolivet if (bnk->pert < bnk->pmin) bnk->pert = 0.0;
652eb910715SAlp Dener break;
653eb910715SAlp Dener }
654fed79b8eSAlp Dener *stepType = BNK_NEWTON;
655eb910715SAlp Dener }
6566b591159SAlp Dener break;
6576b591159SAlp Dener
6586b591159SAlp Dener case BNK_BFGS:
6596b591159SAlp Dener /* Check for success (descent direction) */
6609566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
6616b591159SAlp Dener if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
6626b591159SAlp Dener /* Step is not descent or solve was not successful
6636b591159SAlp Dener Use steepest descent direction (scaled) */
6649566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
6659566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
6669566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection));
6679566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
6689566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
6696b591159SAlp Dener *stepType = BNK_SCALED_GRADIENT;
6706b591159SAlp Dener } else {
6716b591159SAlp Dener *stepType = BNK_BFGS;
6726b591159SAlp Dener }
6736b591159SAlp Dener break;
6746b591159SAlp Dener
675d71ae5a4SJacob Faibussowitsch case BNK_SCALED_GRADIENT:
676d71ae5a4SJacob Faibussowitsch break;
6776b591159SAlp Dener
678d71ae5a4SJacob Faibussowitsch default:
679d71ae5a4SJacob Faibussowitsch break;
6806b591159SAlp Dener }
6813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
682eb910715SAlp Dener }
683eb910715SAlp Dener
684df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
685df278d8fSAlp Dener
686df278d8fSAlp Dener Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
687df278d8fSAlp Dener Newton step does not produce a valid step length.
688df278d8fSAlp Dener */
689df278d8fSAlp Dener
TaoBNKPerformLineSearch(Tao tao,PetscInt * stepType,PetscReal * steplen,TaoLineSearchConvergedReason * reason)690d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
691d71ae5a4SJacob Faibussowitsch {
692c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
693c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason;
694b2d8c577SAlp Dener PetscReal e_min, gdx;
695c14b763aSAlp Dener PetscInt bfgsUpdates;
696c14b763aSAlp Dener
697c14b763aSAlp Dener PetscFunctionBegin;
698c14b763aSAlp Dener /* Perform the linesearch */
6999566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7009566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
701c14b763aSAlp Dener
702b2d8c577SAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
703c14b763aSAlp Dener /* Linesearch failed, revert solution */
704c14b763aSAlp Dener bnk->f = bnk->fold;
7059566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->Xold, tao->solution));
7069566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));
707c14b763aSAlp Dener
708937a31a1SAlp Dener switch (*stepType) {
709c14b763aSAlp Dener case BNK_NEWTON:
7108d5ead36SAlp Dener /* Failed to obtain acceptable iterate with Newton step
711c14b763aSAlp Dener Update the perturbation for next time */
712c14b763aSAlp Dener if (bnk->pert <= 0.0) {
7132e6e4ca1SStefano Zampini PetscBool is_gltr;
7142e6e4ca1SStefano Zampini
715c14b763aSAlp Dener /* Initialize the perturbation */
716c14b763aSAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
717f4f49eeaSPierre Jolivet PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
7182e6e4ca1SStefano Zampini if (is_gltr) {
7199566063dSJacob Faibussowitsch PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
720c14b763aSAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min);
721c14b763aSAlp Dener }
722c14b763aSAlp Dener } else {
723c14b763aSAlp Dener /* Increase the perturbation */
724c14b763aSAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
725c14b763aSAlp Dener }
726c14b763aSAlp Dener
7270ad3a497SAlp Dener if (!bnk->M) {
728c14b763aSAlp Dener /* We don't have the bfgs matrix around and being updated
729c14b763aSAlp Dener Must use gradient direction in this case */
7309566063dSJacob Faibussowitsch PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection));
731937a31a1SAlp Dener *stepType = BNK_GRADIENT;
732c14b763aSAlp Dener } else {
733c14b763aSAlp Dener /* Attempt to use the BFGS direction */
7349566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
7358d5ead36SAlp Dener /* Check for success (descent direction)
7368d5ead36SAlp Dener NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
7379566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
7383105154fSTodd Munson if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
739c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number
740c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case
741c14b763aSAlp Dener Use steepest descent direction (scaled) */
7429566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7439566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7449566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
745c14b763aSAlp Dener
746c14b763aSAlp Dener bfgsUpdates = 1;
747937a31a1SAlp Dener *stepType = BNK_SCALED_GRADIENT;
748c14b763aSAlp Dener } else {
7499566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
750c14b763aSAlp Dener if (1 == bfgsUpdates) {
751c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */
752937a31a1SAlp Dener *stepType = BNK_SCALED_GRADIENT;
753c14b763aSAlp Dener } else {
754937a31a1SAlp Dener *stepType = BNK_BFGS;
755c14b763aSAlp Dener }
756c14b763aSAlp Dener }
757c14b763aSAlp Dener }
758c14b763aSAlp Dener break;
759c14b763aSAlp Dener
760c14b763aSAlp Dener case BNK_BFGS:
761c14b763aSAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS
762c14b763aSAlp Dener Failed to obtain acceptable iterate with BFGS step
763c14b763aSAlp Dener Attempt to use the scaled gradient direction */
7649566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
7659566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
7669566063dSJacob Faibussowitsch PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
767c14b763aSAlp Dener
768c14b763aSAlp Dener bfgsUpdates = 1;
769937a31a1SAlp Dener *stepType = BNK_SCALED_GRADIENT;
770c14b763aSAlp Dener break;
771c14b763aSAlp Dener }
7728d5ead36SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7739566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
7749566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
775c14b763aSAlp Dener
7768d5ead36SAlp Dener /* Perform one last line search with the fall-back step */
7779566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
7789566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
779c14b763aSAlp Dener }
780c14b763aSAlp Dener *reason = ls_reason;
7813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
782c14b763aSAlp Dener }
783c14b763aSAlp Dener
784df278d8fSAlp Dener /* Routine for updating the trust radius.
785df278d8fSAlp Dener
786df278d8fSAlp Dener Function features three different update methods:
787df278d8fSAlp Dener 1) Line-search step length based
788df278d8fSAlp Dener 2) Predicted decrease on the CG quadratic model
789df278d8fSAlp Dener 3) Interpolation
790df278d8fSAlp Dener */
791df278d8fSAlp Dener
TaoBNKUpdateTrustRadius(Tao tao,PetscReal prered,PetscReal actred,PetscInt updateType,PetscInt stepType,PetscBool * accept)792d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
793d71ae5a4SJacob Faibussowitsch {
794080d2917SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
795080d2917SAlp Dener
796b1c2d0e3SAlp Dener PetscReal step, kappa;
797080d2917SAlp Dener PetscReal gdx, tau_1, tau_2, tau_min, tau_max;
798080d2917SAlp Dener
799080d2917SAlp Dener PetscFunctionBegin;
800080d2917SAlp Dener /* Update trust region radius */
801080d2917SAlp Dener *accept = PETSC_FALSE;
80228017e9fSAlp Dener switch (updateType) {
803080d2917SAlp Dener case BNK_UPDATE_STEP:
804c14b763aSAlp Dener *accept = PETSC_TRUE; /* always accept here because line search succeeded */
805080d2917SAlp Dener if (stepType == BNK_NEWTON) {
8069566063dSJacob Faibussowitsch PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step));
807080d2917SAlp Dener if (step < bnk->nu1) {
808080d2917SAlp Dener /* Very bad step taken; reduce radius */
809080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
810080d2917SAlp Dener } else if (step < bnk->nu2) {
811080d2917SAlp Dener /* Reasonably bad step taken; reduce radius */
812080d2917SAlp Dener tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
813080d2917SAlp Dener } else if (step < bnk->nu3) {
814080d2917SAlp Dener /* Reasonable step was taken; leave radius alone */
815080d2917SAlp Dener if (bnk->omega3 < 1.0) {
816080d2917SAlp Dener tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
817080d2917SAlp Dener } else if (bnk->omega3 > 1.0) {
818080d2917SAlp Dener tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
819080d2917SAlp Dener }
820080d2917SAlp Dener } else if (step < bnk->nu4) {
821080d2917SAlp Dener /* Full step taken; increase the radius */
822080d2917SAlp Dener tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
823080d2917SAlp Dener } else {
824080d2917SAlp Dener /* More than full step taken; increase the radius */
825080d2917SAlp Dener tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
826080d2917SAlp Dener }
827080d2917SAlp Dener } else {
828080d2917SAlp Dener /* Newton step was not good; reduce the radius */
829080d2917SAlp Dener tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
830080d2917SAlp Dener }
831080d2917SAlp Dener break;
832080d2917SAlp Dener
833080d2917SAlp Dener case BNK_UPDATE_REDUCTION:
834080d2917SAlp Dener if (stepType == BNK_NEWTON) {
835e0ed867bSAlp Dener if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
836fed79b8eSAlp Dener /* The predicted reduction has the wrong sign. This cannot
837fed79b8eSAlp Dener happen in infinite precision arithmetic. Step should
838fed79b8eSAlp Dener be rejected! */
839080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8403105154fSTodd Munson } else {
841b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) {
842080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
843080d2917SAlp Dener } else {
8443105154fSTodd Munson if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon)) {
845080d2917SAlp Dener kappa = 1.0;
8463105154fSTodd Munson } else {
847080d2917SAlp Dener kappa = actred / prered;
848080d2917SAlp Dener }
849fed79b8eSAlp Dener /* Accept or reject the step and update radius */
850080d2917SAlp Dener if (kappa < bnk->eta1) {
851fed79b8eSAlp Dener /* Reject the step */
852080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
8533105154fSTodd Munson } else {
854fed79b8eSAlp Dener /* Accept the step */
855c133c014SAlp Dener *accept = PETSC_TRUE;
856c133c014SAlp Dener /* Update the trust region radius only if the computed step is at the trust radius boundary */
8578d5ead36SAlp Dener if (bnk->dnorm == tao->trust) {
858080d2917SAlp Dener if (kappa < bnk->eta2) {
859080d2917SAlp Dener /* Marginal bad step */
860c133c014SAlp Dener tao->trust = bnk->alpha2 * tao->trust;
8613105154fSTodd Munson } else if (kappa < bnk->eta3) {
862fed79b8eSAlp Dener /* Reasonable step */
863fed79b8eSAlp Dener tao->trust = bnk->alpha3 * tao->trust;
8643105154fSTodd Munson } else if (kappa < bnk->eta4) {
865080d2917SAlp Dener /* Good step */
866c133c014SAlp Dener tao->trust = bnk->alpha4 * tao->trust;
8673105154fSTodd Munson } else {
868080d2917SAlp Dener /* Very good step */
869c133c014SAlp Dener tao->trust = bnk->alpha5 * tao->trust;
870080d2917SAlp Dener }
871c133c014SAlp Dener }
872080d2917SAlp Dener }
873080d2917SAlp Dener }
874080d2917SAlp Dener }
875080d2917SAlp Dener } else {
876080d2917SAlp Dener /* Newton step was not good; reduce the radius */
877080d2917SAlp Dener tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
878080d2917SAlp Dener }
879080d2917SAlp Dener break;
880080d2917SAlp Dener
881080d2917SAlp Dener default:
882080d2917SAlp Dener if (stepType == BNK_NEWTON) {
883b1c2d0e3SAlp Dener if (prered < 0.0) {
884080d2917SAlp Dener /* The predicted reduction has the wrong sign. This cannot */
885080d2917SAlp Dener /* happen in infinite precision arithmetic. Step should */
886080d2917SAlp Dener /* be rejected! */
887080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
888080d2917SAlp Dener } else {
889b1c2d0e3SAlp Dener if (PetscIsInfOrNanReal(actred)) {
890080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
891080d2917SAlp Dener } else {
892080d2917SAlp Dener if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
893080d2917SAlp Dener kappa = 1.0;
894080d2917SAlp Dener } else {
895080d2917SAlp Dener kappa = actred / prered;
896080d2917SAlp Dener }
897080d2917SAlp Dener
8989566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
899080d2917SAlp Dener tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
900080d2917SAlp Dener tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
901080d2917SAlp Dener tau_min = PetscMin(tau_1, tau_2);
902080d2917SAlp Dener tau_max = PetscMax(tau_1, tau_2);
903080d2917SAlp Dener
904080d2917SAlp Dener if (kappa >= 1.0 - bnk->mu1) {
905080d2917SAlp Dener /* Great agreement */
906080d2917SAlp Dener *accept = PETSC_TRUE;
907080d2917SAlp Dener if (tau_max < 1.0) {
908080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
909080d2917SAlp Dener } else if (tau_max > bnk->gamma4) {
910080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
911080d2917SAlp Dener } else {
912080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
913080d2917SAlp Dener }
914080d2917SAlp Dener } else if (kappa >= 1.0 - bnk->mu2) {
915080d2917SAlp Dener /* Good agreement */
916080d2917SAlp Dener *accept = PETSC_TRUE;
917080d2917SAlp Dener if (tau_max < bnk->gamma2) {
918080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
919080d2917SAlp Dener } else if (tau_max > bnk->gamma3) {
920080d2917SAlp Dener tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
921080d2917SAlp Dener } else if (tau_max < 1.0) {
922080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
923080d2917SAlp Dener } else {
924080d2917SAlp Dener tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
925080d2917SAlp Dener }
926080d2917SAlp Dener } else {
927080d2917SAlp Dener /* Not good agreement */
928080d2917SAlp Dener if (tau_min > 1.0) {
929080d2917SAlp Dener tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
930080d2917SAlp Dener } else if (tau_max < bnk->gamma1) {
931080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
932080d2917SAlp Dener } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
933080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
934080d2917SAlp Dener } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
935080d2917SAlp Dener tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
936080d2917SAlp Dener } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
937080d2917SAlp Dener tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
938080d2917SAlp Dener } else {
939080d2917SAlp Dener tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
940080d2917SAlp Dener }
941080d2917SAlp Dener }
942080d2917SAlp Dener }
943080d2917SAlp Dener }
944080d2917SAlp Dener } else {
945080d2917SAlp Dener /* Newton step was not good; reduce the radius */
946080d2917SAlp Dener tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
947080d2917SAlp Dener }
94828017e9fSAlp Dener break;
949080d2917SAlp Dener }
950c133c014SAlp Dener /* Make sure the radius does not violate min and max settings */
951c133c014SAlp Dener tao->trust = PetscMin(tao->trust, bnk->max_radius);
952fed79b8eSAlp Dener tao->trust = PetscMax(tao->trust, bnk->min_radius);
9533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
954080d2917SAlp Dener }
955080d2917SAlp Dener
TaoBNKAddStepCounts(Tao tao,PetscInt stepType)956d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
957d71ae5a4SJacob Faibussowitsch {
95862675beeSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
95962675beeSAlp Dener
96062675beeSAlp Dener PetscFunctionBegin;
96162675beeSAlp Dener switch (stepType) {
962d71ae5a4SJacob Faibussowitsch case BNK_NEWTON:
963d71ae5a4SJacob Faibussowitsch ++bnk->newt;
964d71ae5a4SJacob Faibussowitsch break;
965d71ae5a4SJacob Faibussowitsch case BNK_BFGS:
966d71ae5a4SJacob Faibussowitsch ++bnk->bfgs;
967d71ae5a4SJacob Faibussowitsch break;
968d71ae5a4SJacob Faibussowitsch case BNK_SCALED_GRADIENT:
969d71ae5a4SJacob Faibussowitsch ++bnk->sgrad;
970d71ae5a4SJacob Faibussowitsch break;
971d71ae5a4SJacob Faibussowitsch case BNK_GRADIENT:
972d71ae5a4SJacob Faibussowitsch ++bnk->grad;
973d71ae5a4SJacob Faibussowitsch break;
974d71ae5a4SJacob Faibussowitsch default:
975d71ae5a4SJacob Faibussowitsch break;
97662675beeSAlp Dener }
9773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97862675beeSAlp Dener }
97962675beeSAlp Dener
TaoSetUp_BNK(Tao tao)980d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetUp_BNK(Tao tao)
981d71ae5a4SJacob Faibussowitsch {
982eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
983eb910715SAlp Dener
984eb910715SAlp Dener PetscFunctionBegin;
98548a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
98648a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
98748a46eb9SPierre Jolivet if (!bnk->W) PetscCall(VecDuplicate(tao->solution, &bnk->W));
98848a46eb9SPierre Jolivet if (!bnk->Xold) PetscCall(VecDuplicate(tao->solution, &bnk->Xold));
98948a46eb9SPierre Jolivet if (!bnk->Gold) PetscCall(VecDuplicate(tao->solution, &bnk->Gold));
99048a46eb9SPierre Jolivet if (!bnk->Xwork) PetscCall(VecDuplicate(tao->solution, &bnk->Xwork));
99148a46eb9SPierre Jolivet if (!bnk->Gwork) PetscCall(VecDuplicate(tao->solution, &bnk->Gwork));
99248a46eb9SPierre Jolivet if (!bnk->unprojected_gradient) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient));
99348a46eb9SPierre Jolivet if (!bnk->unprojected_gradient_old) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient_old));
99448a46eb9SPierre Jolivet if (!bnk->Diag_min) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_min));
99548a46eb9SPierre Jolivet if (!bnk->Diag_max) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_max));
996e031d6f5SAlp Dener if (bnk->max_cg_its > 0) {
997c4b75bccSAlp Dener /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
998c4b75bccSAlp Dener bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
999f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)bnk->unprojected_gradient_old));
10009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old));
100189da521bSAlp Dener bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1002f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)bnk->unprojected_gradient));
10039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient));
1004c4b75bccSAlp Dener bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1005f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)bnk->Gold));
10069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg_ctx->G_old));
1007c4b75bccSAlp Dener bnk->bncg_ctx->G_old = bnk->Gold;
1008f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)tao->gradient));
10099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg->gradient));
1010c4b75bccSAlp Dener bnk->bncg->gradient = tao->gradient;
1011f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)tao->stepdirection));
10129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->bncg->stepdirection));
1013c4b75bccSAlp Dener bnk->bncg->stepdirection = tao->stepdirection;
10149566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(bnk->bncg, tao->solution));
1015c4b75bccSAlp Dener /* Copy over some settings from BNK into BNCG */
10169566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its));
10179566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol));
10189566063dSJacob Faibussowitsch PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin));
10199566063dSJacob Faibussowitsch PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP));
10209566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(bnk->bncg, tao->ops->computeobjective, tao->user_objP));
10219566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(bnk->bncg, NULL, tao->ops->computegradient, tao->user_gradP));
10229566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(bnk->bncg, NULL, tao->ops->computeobjectiveandgradient, tao->user_objgradP));
1023f4f49eeaSPierre Jolivet PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)bnk->bncg));
1024e031d6f5SAlp Dener }
102583c8fe1dSLisandro Dalcin bnk->X_inactive = NULL;
102683c8fe1dSLisandro Dalcin bnk->G_inactive = NULL;
102783c8fe1dSLisandro Dalcin bnk->inactive_work = NULL;
102883c8fe1dSLisandro Dalcin bnk->active_work = NULL;
102983c8fe1dSLisandro Dalcin bnk->inactive_idx = NULL;
103083c8fe1dSLisandro Dalcin bnk->active_idx = NULL;
103183c8fe1dSLisandro Dalcin bnk->active_lower = NULL;
103283c8fe1dSLisandro Dalcin bnk->active_upper = NULL;
103383c8fe1dSLisandro Dalcin bnk->active_fixed = NULL;
103483c8fe1dSLisandro Dalcin bnk->M = NULL;
103583c8fe1dSLisandro Dalcin bnk->H_inactive = NULL;
103683c8fe1dSLisandro Dalcin bnk->Hpre_inactive = NULL;
10373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1038eb910715SAlp Dener }
1039eb910715SAlp Dener
TaoDestroy_BNK(Tao tao)1040d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDestroy_BNK(Tao tao)
1041d71ae5a4SJacob Faibussowitsch {
1042eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
1043eb910715SAlp Dener
1044eb910715SAlp Dener PetscFunctionBegin;
10459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->W));
10469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Xold));
10479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Gold));
10489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Xwork));
10499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Gwork));
10509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->unprojected_gradient));
10519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->unprojected_gradient_old));
10529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Diag_min));
10539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bnk->Diag_max));
10549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_lower));
10559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_upper));
10569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_fixed));
10579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->active_idx));
10589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bnk->inactive_idx));
10599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->Hpre_inactive));
10609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&bnk->H_inactive));
10619566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&bnk->bncg));
1062a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp));
10639566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
10643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1065eb910715SAlp Dener }
1066eb910715SAlp Dener
TaoSetFromOptions_BNK(Tao tao,PetscOptionItems PetscOptionsObject)1067ce78bad3SBarry Smith PetscErrorCode TaoSetFromOptions_BNK(Tao tao, PetscOptionItems PetscOptionsObject)
1068d71ae5a4SJacob Faibussowitsch {
1069eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
1070eb910715SAlp Dener
1071eb910715SAlp Dener PetscFunctionBegin;
1072d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Newton-Krylov method for bound constrained optimization");
10739566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL));
10749566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL));
10759566063dSJacob 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));
10769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval, NULL));
10779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin, NULL));
10789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax, NULL));
10799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac, NULL));
10809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin, NULL));
10819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax, NULL));
10829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac, NULL));
10839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac, NULL));
10849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac, NULL));
10859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac, NULL));
10869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1, NULL));
10879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2, NULL));
10889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3, NULL));
10899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4, NULL));
10909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1, NULL));
10919566063dSJacob 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));
10929566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3, NULL));
10939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4, NULL));
10949566063dSJacob 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));
10959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1, NULL));
10969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2, NULL));
10979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3, NULL));
10989566063dSJacob 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));
10999566063dSJacob 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));
11009566063dSJacob 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));
11019566063dSJacob 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));
11029566063dSJacob 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));
11039566063dSJacob 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));
11049566063dSJacob 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));
11059566063dSJacob 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));
11069566063dSJacob 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));
11079566063dSJacob 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));
11089566063dSJacob 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));
11099566063dSJacob 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));
11109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i, NULL));
11119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1, NULL));
11129566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2, NULL));
11139566063dSJacob 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));
11149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2, NULL));
11159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3, NULL));
11169566063dSJacob 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));
11179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta, NULL));
11189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius, NULL));
11199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius, NULL));
11209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
11219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
11229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
11239566063dSJacob 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));
1124d0609cedSBarry Smith PetscOptionsHeadEnd();
11258ebe3e4eSStefano Zampini
1126f4f49eeaSPierre Jolivet PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)tao)->prefix));
11279566063dSJacob Faibussowitsch PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_cg_"));
11289566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(bnk->bncg));
11298ebe3e4eSStefano Zampini
1130f4f49eeaSPierre Jolivet PetscCall(KSPSetOptionsPrefix(tao->ksp, ((PetscObject)tao)->prefix));
11319566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_bnk_"));
11329566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp));
11333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1134eb910715SAlp Dener }
1135eb910715SAlp Dener
TaoView_BNK(Tao tao,PetscViewer viewer)1136d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1137d71ae5a4SJacob Faibussowitsch {
1138eb910715SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
1139eb910715SAlp Dener PetscInt nrejects;
1140eb910715SAlp Dener PetscBool isascii;
1141eb910715SAlp Dener
1142eb910715SAlp Dener PetscFunctionBegin;
11439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1144eb910715SAlp Dener if (isascii) {
11459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
1146b17ffb64SBarry Smith PetscCall(TaoView(bnk->bncg, viewer));
1147b9ac7092SAlp Dener if (bnk->M) {
11489566063dSJacob Faibussowitsch PetscCall(MatLMVMGetRejectCount(bnk->M, &nrejects));
114963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n", nrejects));
1150eb910715SAlp Dener }
115163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its));
115263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt));
115348a46eb9SPierre Jolivet if (bnk->M) PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs));
115463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad));
115563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad));
11569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n"));
115763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " atol: %" PetscInt_FMT "\n", bnk->ksp_atol));
115863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol));
115963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol));
116063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " negc: %" PetscInt_FMT "\n", bnk->ksp_negc));
116163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol));
116263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " iter: %" PetscInt_FMT "\n", bnk->ksp_iter));
116363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " othr: %" PetscInt_FMT "\n", bnk->ksp_othr));
11649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
1165eb910715SAlp Dener }
11663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1167eb910715SAlp Dener }
1168eb910715SAlp Dener
1169eb910715SAlp Dener /*MC
1170eb910715SAlp Dener TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
117166ed3702SAlp Dener At each iteration, the BNK methods solve the symmetric
11721cb3f120SPierre Jolivet system of equations to obtain the step direction dk:
1173eb910715SAlp Dener Hk dk = -gk
11742b97c8d8SAlp Dener for free variables only. The step can be globalized either through
11752b97c8d8SAlp Dener trust-region methods, or a line search, or a heuristic mixture of both.
1176eb910715SAlp Dener
1177eb910715SAlp Dener Options Database Keys:
11789fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
11799fa2b5dcSStefano Zampini . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
11809fa2b5dcSStefano Zampini . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
11819fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
11829fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
11839fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
11849fa2b5dcSStefano Zampini . -tao_bnk_sval - (developer) Hessian perturbation starting value
11859fa2b5dcSStefano Zampini . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
11869fa2b5dcSStefano Zampini . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
11879fa2b5dcSStefano Zampini . -tao_bnk_pmin - (developer) minimum Hessian perturbation
11889fa2b5dcSStefano Zampini . -tao_bnk_pmax - (developer) aximum Hessian perturbation
11899fa2b5dcSStefano Zampini . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
11909fa2b5dcSStefano Zampini . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
11919fa2b5dcSStefano Zampini . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
11929fa2b5dcSStefano Zampini . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
11939fa2b5dcSStefano Zampini . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
11949fa2b5dcSStefano Zampini . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction)
11959fa2b5dcSStefano Zampini . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
11969fa2b5dcSStefano Zampini . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
11979fa2b5dcSStefano Zampini . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction)
11989fa2b5dcSStefano Zampini . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
11999fa2b5dcSStefano Zampini . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
12009fa2b5dcSStefano Zampini . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
12019fa2b5dcSStefano Zampini . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
12029fa2b5dcSStefano Zampini . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
12039fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
12049fa2b5dcSStefano Zampini . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
12059fa2b5dcSStefano Zampini . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation)
12069fa2b5dcSStefano Zampini . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
12079fa2b5dcSStefano Zampini . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
12089fa2b5dcSStefano Zampini . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
12099fa2b5dcSStefano Zampini . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
12109fa2b5dcSStefano Zampini . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation)
12119fa2b5dcSStefano Zampini . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step)
12129fa2b5dcSStefano Zampini . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
12139fa2b5dcSStefano Zampini . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step)
12149fa2b5dcSStefano Zampini . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step)
12159fa2b5dcSStefano Zampini . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
12169fa2b5dcSStefano Zampini . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
12179fa2b5dcSStefano Zampini . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step)
12189fa2b5dcSStefano Zampini . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
12199fa2b5dcSStefano Zampini . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
12209fa2b5dcSStefano Zampini . -tao_bnk_mu1_i - (developer) threshold for accepting very good step (-init_type interpolation)
12219fa2b5dcSStefano Zampini . -tao_bnk_mu2_i - (developer) threshold for accepting good step (-init_type interpolation)
12229fa2b5dcSStefano Zampini . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
12239fa2b5dcSStefano Zampini . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
12249fa2b5dcSStefano Zampini . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
12259fa2b5dcSStefano Zampini . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
12269fa2b5dcSStefano Zampini - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation)
1227eb910715SAlp Dener
1228eb910715SAlp Dener Level: beginner
1229eb910715SAlp Dener M*/
1230eb910715SAlp Dener
TaoCreate_BNK(Tao tao)1231d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoCreate_BNK(Tao tao)
1232d71ae5a4SJacob Faibussowitsch {
1233eb910715SAlp Dener TAO_BNK *bnk;
1234b9ac7092SAlp Dener PC pc;
1235eb910715SAlp Dener
1236eb910715SAlp Dener PetscFunctionBegin;
12374dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bnk));
1238eb910715SAlp Dener
1239eb910715SAlp Dener tao->ops->setup = TaoSetUp_BNK;
1240eb910715SAlp Dener tao->ops->view = TaoView_BNK;
1241eb910715SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1242eb910715SAlp Dener tao->ops->destroy = TaoDestroy_BNK;
1243eb910715SAlp Dener
1244eb910715SAlp Dener /* Override default settings (unless already changed) */
1245606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
1246606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 50);
1247606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, trust0, 100.0);
1248eb910715SAlp Dener
1249eb910715SAlp Dener tao->data = (void *)bnk;
1250eb910715SAlp Dener
125166ed3702SAlp Dener /* Hessian shifting parameters */
1252e0ed867bSAlp Dener bnk->computehessian = TaoBNKComputeHessian;
1253e0ed867bSAlp Dener bnk->computestep = TaoBNKComputeStep;
1254e0ed867bSAlp Dener
1255eb910715SAlp Dener bnk->sval = 0.0;
1256eb910715SAlp Dener bnk->imin = 1.0e-4;
1257eb910715SAlp Dener bnk->imax = 1.0e+2;
1258eb910715SAlp Dener bnk->imfac = 1.0e-1;
1259eb910715SAlp Dener
1260eb910715SAlp Dener bnk->pmin = 1.0e-12;
1261eb910715SAlp Dener bnk->pmax = 1.0e+2;
1262eb910715SAlp Dener bnk->pgfac = 1.0e+1;
1263eb910715SAlp Dener bnk->psfac = 4.0e-1;
1264eb910715SAlp Dener bnk->pmgfac = 1.0e-1;
1265eb910715SAlp Dener bnk->pmsfac = 1.0e-1;
1266eb910715SAlp Dener
1267eb910715SAlp Dener /* Default values for trust-region radius update based on steplength */
1268eb910715SAlp Dener bnk->nu1 = 0.25;
1269eb910715SAlp Dener bnk->nu2 = 0.50;
1270eb910715SAlp Dener bnk->nu3 = 1.00;
1271eb910715SAlp Dener bnk->nu4 = 1.25;
1272eb910715SAlp Dener
1273eb910715SAlp Dener bnk->omega1 = 0.25;
1274eb910715SAlp Dener bnk->omega2 = 0.50;
1275eb910715SAlp Dener bnk->omega3 = 1.00;
1276eb910715SAlp Dener bnk->omega4 = 2.00;
1277eb910715SAlp Dener bnk->omega5 = 4.00;
1278eb910715SAlp Dener
1279eb910715SAlp Dener /* Default values for trust-region radius update based on reduction */
1280eb910715SAlp Dener bnk->eta1 = 1.0e-4;
1281eb910715SAlp Dener bnk->eta2 = 0.25;
1282eb910715SAlp Dener bnk->eta3 = 0.50;
1283eb910715SAlp Dener bnk->eta4 = 0.90;
1284eb910715SAlp Dener
1285eb910715SAlp Dener bnk->alpha1 = 0.25;
1286eb910715SAlp Dener bnk->alpha2 = 0.50;
1287eb910715SAlp Dener bnk->alpha3 = 1.00;
1288eb910715SAlp Dener bnk->alpha4 = 2.00;
1289eb910715SAlp Dener bnk->alpha5 = 4.00;
1290eb910715SAlp Dener
1291eb910715SAlp Dener /* Default values for trust-region radius update based on interpolation */
1292eb910715SAlp Dener bnk->mu1 = 0.10;
1293eb910715SAlp Dener bnk->mu2 = 0.50;
1294eb910715SAlp Dener
1295eb910715SAlp Dener bnk->gamma1 = 0.25;
1296eb910715SAlp Dener bnk->gamma2 = 0.50;
1297eb910715SAlp Dener bnk->gamma3 = 2.00;
1298eb910715SAlp Dener bnk->gamma4 = 4.00;
1299eb910715SAlp Dener
1300eb910715SAlp Dener bnk->theta = 0.05;
1301eb910715SAlp Dener
1302eb910715SAlp Dener /* Default values for trust region initialization based on interpolation */
1303eb910715SAlp Dener bnk->mu1_i = 0.35;
1304eb910715SAlp Dener bnk->mu2_i = 0.50;
1305eb910715SAlp Dener
1306eb910715SAlp Dener bnk->gamma1_i = 0.0625;
1307eb910715SAlp Dener bnk->gamma2_i = 0.5;
1308eb910715SAlp Dener bnk->gamma3_i = 2.0;
1309eb910715SAlp Dener bnk->gamma4_i = 5.0;
1310eb910715SAlp Dener
1311eb910715SAlp Dener bnk->theta_i = 0.25;
1312eb910715SAlp Dener
1313eb910715SAlp Dener /* Remaining parameters */
1314c0f10754SAlp Dener bnk->max_cg_its = 0;
1315eb910715SAlp Dener bnk->min_radius = 1.0e-10;
1316eb910715SAlp Dener bnk->max_radius = 1.0e10;
1317770b7498SAlp Dener bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0);
13180a4511e9SAlp Dener bnk->as_tol = 1.0e-3;
13190a4511e9SAlp Dener bnk->as_step = 1.0e-3;
132062675beeSAlp Dener bnk->dmin = 1.0e-6;
132162675beeSAlp Dener bnk->dmax = 1.0e6;
1322eb910715SAlp Dener
132383c8fe1dSLisandro Dalcin bnk->M = NULL;
132483c8fe1dSLisandro Dalcin bnk->bfgs_pre = NULL;
1325eb910715SAlp Dener bnk->init_type = BNK_INIT_INTERPOLATION;
13267b1c7716SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION;
13272f75a4aaSAlp Dener bnk->as_type = BNK_AS_BERTSEKAS;
1328eb910715SAlp Dener
1329e031d6f5SAlp Dener /* Create the embedded BNCG solver */
13309566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg));
13319566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1));
13329566063dSJacob Faibussowitsch PetscCall(TaoSetType(bnk->bncg, TAOBNCG));
1333e031d6f5SAlp Dener
1334c0f10754SAlp Dener /* Create the line search */
13359566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
13369566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1337f4db9bf7SStefano Zampini PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
13389566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
1339eb910715SAlp Dener
1340eb910715SAlp Dener /* Set linear solver to default for symmetric matrices */
13419566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
13429566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
13439566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp, KSPSTCG));
13449566063dSJacob Faibussowitsch PetscCall(KSPGetPC(tao->ksp, &pc));
13459566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCLMVM));
13463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1347eb910715SAlp Dener }
1348