1 #include <petsctaolinesearch.h> 2 #include <../src/tao/bound/impls/bnk/bnk.h> 3 #include <petscksp.h> 4 5 static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"}; 6 static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"}; 7 static const char *BNK_AS[64] = {"none", "bertsekas"}; 8 9 /*------------------------------------------------------------*/ 10 11 /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */ 12 13 PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH) { 14 TAO_BNK *bnk = (TAO_BNK *)tao->data; 15 PC pc; 16 PetscReal f_min, ftrial, prered, actred, kappa, sigma, resnorm; 17 PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 18 PetscBool is_bfgs, is_jacobi, is_symmetric, sym_set; 19 PetscInt n, N, nDiff; 20 PetscInt i_max = 5; 21 PetscInt j_max = 1; 22 PetscInt i, j; 23 PetscVoidFunction kspTR; 24 25 PetscFunctionBegin; 26 /* Project the current point onto the feasible set */ 27 PetscCall(TaoComputeVariableBounds(tao)); 28 PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU)); 29 if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU)); 30 31 /* Project the initial point onto the feasible region */ 32 PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution)); 33 34 /* Check convergence criteria */ 35 PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient)); 36 PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 37 PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient)); 38 PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0)); 39 PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm)); 40 41 /* Test the initial point for convergence */ 42 PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W)); 43 PetscCall(VecNorm(bnk->W, NORM_2, &resnorm)); 44 PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 45 PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its)); 46 PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0)); 47 PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 48 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 49 50 /* Reset KSP stopping reason counters */ 51 bnk->ksp_atol = 0; 52 bnk->ksp_rtol = 0; 53 bnk->ksp_dtol = 0; 54 bnk->ksp_ctol = 0; 55 bnk->ksp_negc = 0; 56 bnk->ksp_iter = 0; 57 bnk->ksp_othr = 0; 58 59 /* Reset accepted step type counters */ 60 bnk->tot_cg_its = 0; 61 bnk->newt = 0; 62 bnk->bfgs = 0; 63 bnk->sgrad = 0; 64 bnk->grad = 0; 65 66 /* Initialize the Hessian perturbation */ 67 bnk->pert = bnk->sval; 68 69 /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */ 70 PetscCall(VecSet(tao->stepdirection, 0.0)); 71 72 /* Allocate the vectors needed for the BFGS approximation */ 73 PetscCall(KSPGetPC(tao->ksp, &pc)); 74 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs)); 75 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi)); 76 if (is_bfgs) { 77 bnk->bfgs_pre = pc; 78 PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M)); 79 PetscCall(VecGetLocalSize(tao->solution, &n)); 80 PetscCall(VecGetSize(tao->solution, &N)); 81 PetscCall(MatSetSizes(bnk->M, n, n, N, N)); 82 PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient)); 83 PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric)); 84 PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric."); 85 } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE)); 86 87 /* Prepare the min/max vectors for safeguarding diagonal scales */ 88 PetscCall(VecSet(bnk->Diag_min, bnk->dmin)); 89 PetscCall(VecSet(bnk->Diag_max, bnk->dmax)); 90 91 /* Initialize trust-region radius. The initialization is only performed 92 when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 93 *needH = PETSC_TRUE; 94 PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp, "KSPCGSetRadius_C", &kspTR)); 95 if (kspTR) { 96 switch (initType) { 97 case BNK_INIT_CONSTANT: 98 /* Use the initial radius specified */ 99 tao->trust = tao->trust0; 100 break; 101 102 case BNK_INIT_INTERPOLATION: 103 /* Use interpolation based on the initial Hessian */ 104 max_radius = 0.0; 105 tao->trust = tao->trust0; 106 for (j = 0; j < j_max; ++j) { 107 f_min = bnk->f; 108 sigma = 0.0; 109 110 if (*needH) { 111 /* Compute the Hessian at the new step, and extract the inactive subsystem */ 112 PetscCall((*bnk->computehessian)(tao)); 113 PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE)); 114 PetscCall(MatDestroy(&bnk->H_inactive)); 115 if (bnk->active_idx) { 116 PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive)); 117 } else { 118 PetscCall(PetscObjectReference((PetscObject)tao->hessian)); 119 bnk->H_inactive = tao->hessian; 120 } 121 *needH = PETSC_FALSE; 122 } 123 124 for (i = 0; i < i_max; ++i) { 125 /* Take a steepest descent step and snap it to bounds */ 126 PetscCall(VecCopy(tao->solution, bnk->Xold)); 127 PetscCall(VecAXPY(tao->solution, -tao->trust / bnk->gnorm, tao->gradient)); 128 PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution)); 129 /* Compute the step we actually accepted */ 130 PetscCall(VecCopy(tao->solution, bnk->W)); 131 PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold)); 132 /* Compute the objective at the trial */ 133 PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial)); 134 PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 135 PetscCall(VecCopy(bnk->Xold, tao->solution)); 136 if (PetscIsInfOrNanReal(ftrial)) { 137 tau = bnk->gamma1_i; 138 } else { 139 if (ftrial < f_min) { 140 f_min = ftrial; 141 sigma = -tao->trust / bnk->gnorm; 142 } 143 144 /* Compute the predicted and actual reduction */ 145 if (bnk->active_idx) { 146 PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive)); 147 PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 148 } else { 149 bnk->X_inactive = bnk->W; 150 bnk->inactive_work = bnk->Xwork; 151 } 152 PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work)); 153 PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered)); 154 if (bnk->active_idx) { 155 PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive)); 156 PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 157 } 158 prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 159 actred = bnk->f - ftrial; 160 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 161 kappa = 1.0; 162 } else { 163 kappa = actred / prered; 164 } 165 166 tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 167 tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 168 tau_min = PetscMin(tau_1, tau_2); 169 tau_max = PetscMax(tau_1, tau_2); 170 171 if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) { 172 /* Great agreement */ 173 max_radius = PetscMax(max_radius, tao->trust); 174 175 if (tau_max < 1.0) { 176 tau = bnk->gamma3_i; 177 } else if (tau_max > bnk->gamma4_i) { 178 tau = bnk->gamma4_i; 179 } else { 180 tau = tau_max; 181 } 182 } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) { 183 /* Good agreement */ 184 max_radius = PetscMax(max_radius, tao->trust); 185 186 if (tau_max < bnk->gamma2_i) { 187 tau = bnk->gamma2_i; 188 } else if (tau_max > bnk->gamma3_i) { 189 tau = bnk->gamma3_i; 190 } else { 191 tau = tau_max; 192 } 193 } else { 194 /* Not good agreement */ 195 if (tau_min > 1.0) { 196 tau = bnk->gamma2_i; 197 } else if (tau_max < bnk->gamma1_i) { 198 tau = bnk->gamma1_i; 199 } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 200 tau = bnk->gamma1_i; 201 } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 202 tau = tau_1; 203 } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 204 tau = tau_2; 205 } else { 206 tau = tau_max; 207 } 208 } 209 } 210 tao->trust = tau * tao->trust; 211 } 212 213 if (f_min < bnk->f) { 214 /* We accidentally found a solution better than the initial, so accept it */ 215 bnk->f = f_min; 216 PetscCall(VecCopy(tao->solution, bnk->Xold)); 217 PetscCall(VecAXPY(tao->solution, sigma, tao->gradient)); 218 PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution)); 219 PetscCall(VecCopy(tao->solution, tao->stepdirection)); 220 PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold)); 221 PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient)); 222 PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 223 PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient)); 224 PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0)); 225 /* Compute gradient at the new iterate and flip switch to compute the Hessian later */ 226 PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm)); 227 *needH = PETSC_TRUE; 228 /* Test the new step for convergence */ 229 PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W)); 230 PetscCall(VecNorm(bnk->W, NORM_2, &resnorm)); 231 PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 232 PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its)); 233 PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0)); 234 PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 235 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 236 /* active BNCG recycling early because we have a stepdirection computed */ 237 PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE)); 238 } 239 } 240 tao->trust = PetscMax(tao->trust, max_radius); 241 242 /* Ensure that the trust radius is within the limits */ 243 tao->trust = PetscMax(tao->trust, bnk->min_radius); 244 tao->trust = PetscMin(tao->trust, bnk->max_radius); 245 break; 246 247 default: 248 /* Norm of the first direction will initialize radius */ 249 tao->trust = 0.0; 250 break; 251 } 252 } 253 PetscFunctionReturn(0); 254 } 255 256 /*------------------------------------------------------------*/ 257 258 /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */ 259 260 PetscErrorCode TaoBNKComputeHessian(Tao tao) { 261 TAO_BNK *bnk = (TAO_BNK *)tao->data; 262 263 PetscFunctionBegin; 264 /* Compute the Hessian */ 265 PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre)); 266 /* Add a correction to the BFGS preconditioner */ 267 if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 268 /* Prepare the reduced sub-matrices for the inactive set */ 269 PetscCall(MatDestroy(&bnk->Hpre_inactive)); 270 PetscCall(MatDestroy(&bnk->H_inactive)); 271 if (bnk->active_idx) { 272 PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive)); 273 if (tao->hessian == tao->hessian_pre) { 274 PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive)); 275 bnk->Hpre_inactive = bnk->H_inactive; 276 } else { 277 PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive)); 278 } 279 if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx)); 280 } else { 281 PetscCall(PetscObjectReference((PetscObject)tao->hessian)); 282 bnk->H_inactive = tao->hessian; 283 if (tao->hessian == tao->hessian_pre) { 284 PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive)); 285 bnk->Hpre_inactive = bnk->H_inactive; 286 } else { 287 PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre)); 288 bnk->Hpre_inactive = tao->hessian_pre; 289 } 290 if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre)); 291 } 292 PetscFunctionReturn(0); 293 } 294 295 /*------------------------------------------------------------*/ 296 297 /* Routine for estimating the active set */ 298 299 PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType) { 300 TAO_BNK *bnk = (TAO_BNK *)tao->data; 301 PetscBool hessComputed, diagExists, hadactive; 302 303 PetscFunctionBegin; 304 hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE; 305 switch (asType) { 306 case BNK_AS_NONE: 307 PetscCall(ISDestroy(&bnk->inactive_idx)); 308 PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx)); 309 PetscCall(ISDestroy(&bnk->active_idx)); 310 PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx)); 311 break; 312 313 case BNK_AS_BERTSEKAS: 314 /* Compute the trial step vector with which we will estimate the active set at the next iteration */ 315 if (bnk->M) { 316 /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */ 317 PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W)); 318 } else { 319 hessComputed = diagExists = PETSC_FALSE; 320 if (tao->hessian) { PetscCall(MatAssembled(tao->hessian, &hessComputed)); } 321 if (hessComputed) { PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists)); } 322 if (diagExists) { 323 /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */ 324 PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork)); 325 PetscCall(VecAbs(bnk->Xwork)); 326 PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork)); 327 PetscCall(VecReciprocal(bnk->Xwork)); 328 PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient)); 329 } else { 330 /* If the Hessian or its diagonal does not exist, we will simply use gradient step */ 331 PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W)); 332 } 333 } 334 PetscCall(VecScale(bnk->W, -1.0)); 335 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)); 336 break; 337 338 default: break; 339 } 340 bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */ 341 PetscFunctionReturn(0); 342 } 343 344 /*------------------------------------------------------------*/ 345 346 /* Routine for bounding the step direction */ 347 348 PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step) { 349 TAO_BNK *bnk = (TAO_BNK *)tao->data; 350 351 PetscFunctionBegin; 352 switch (asType) { 353 case BNK_AS_NONE: PetscCall(VecISSet(step, bnk->active_idx, 0.0)); break; 354 355 case BNK_AS_BERTSEKAS: PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step)); break; 356 357 default: break; 358 } 359 PetscFunctionReturn(0); 360 } 361 362 /*------------------------------------------------------------*/ 363 364 /* Routine for taking a finite number of BNCG iterations to 365 accelerate Newton convergence. 366 367 In practice, this approach simply trades off Hessian evaluations 368 for more gradient evaluations. 369 */ 370 371 PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate) { 372 TAO_BNK *bnk = (TAO_BNK *)tao->data; 373 374 PetscFunctionBegin; 375 *terminate = PETSC_FALSE; 376 if (bnk->max_cg_its > 0) { 377 /* Copy the current function value (important vectors are already shared) */ 378 bnk->bncg_ctx->f = bnk->f; 379 /* Take some small finite number of BNCG iterations */ 380 PetscCall(TaoSolve(bnk->bncg)); 381 /* Add the number of gradient and function evaluations to the total */ 382 tao->nfuncs += bnk->bncg->nfuncs; 383 tao->nfuncgrads += bnk->bncg->nfuncgrads; 384 tao->ngrads += bnk->bncg->ngrads; 385 tao->nhess += bnk->bncg->nhess; 386 bnk->tot_cg_its += bnk->bncg->niter; 387 /* Extract the BNCG function value out and save it into BNK */ 388 bnk->f = bnk->bncg_ctx->f; 389 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) { 390 *terminate = PETSC_TRUE; 391 } else { 392 PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type)); 393 } 394 } 395 PetscFunctionReturn(0); 396 } 397 398 /*------------------------------------------------------------*/ 399 400 /* Routine for computing the Newton step. */ 401 402 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) { 403 TAO_BNK *bnk = (TAO_BNK *)tao->data; 404 PetscInt bfgsUpdates = 0; 405 PetscInt kspits; 406 PetscBool is_lmvm; 407 PetscVoidFunction kspTR; 408 409 PetscFunctionBegin; 410 /* If there are no inactive variables left, save some computation and return an adjusted zero step 411 that has (l-x) and (u-x) for lower and upper bounded variables. */ 412 if (!bnk->inactive_idx) { 413 PetscCall(VecSet(tao->stepdirection, 0.0)); 414 PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 415 PetscFunctionReturn(0); 416 } 417 418 /* Shift the reduced Hessian matrix */ 419 if (shift && bnk->pert > 0) { 420 PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm)); 421 if (is_lmvm) { 422 PetscCall(MatShift(tao->hessian, bnk->pert)); 423 } else { 424 PetscCall(MatShift(bnk->H_inactive, bnk->pert)); 425 if (bnk->H_inactive != bnk->Hpre_inactive) { PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert)); } 426 } 427 } 428 429 /* Solve the Newton system of equations */ 430 tao->ksp_its = 0; 431 PetscCall(VecSet(tao->stepdirection, 0.0)); 432 if (bnk->resetksp) { 433 PetscCall(KSPReset(tao->ksp)); 434 PetscCall(KSPResetFromOptions(tao->ksp)); 435 bnk->resetksp = PETSC_FALSE; 436 } 437 PetscCall(KSPSetOperators(tao->ksp, bnk->H_inactive, bnk->Hpre_inactive)); 438 PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork)); 439 if (bnk->active_idx) { 440 PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 441 PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 442 } else { 443 bnk->G_inactive = bnk->unprojected_gradient; 444 bnk->X_inactive = tao->stepdirection; 445 } 446 PetscCall(KSPCGSetRadius(tao->ksp, tao->trust)); 447 PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive)); 448 PetscCall(KSPGetIterationNumber(tao->ksp, &kspits)); 449 tao->ksp_its += kspits; 450 tao->ksp_tot_its += kspits; 451 PetscCall(PetscObjectQueryFunction((PetscObject)tao->ksp, "KSPCGGetNormD_C", &kspTR)); 452 if (kspTR) { 453 PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm)); 454 455 if (0.0 == tao->trust) { 456 /* Radius was uninitialized; use the norm of the direction */ 457 if (bnk->dnorm > 0.0) { 458 tao->trust = bnk->dnorm; 459 460 /* Modify the radius if it is too large or small */ 461 tao->trust = PetscMax(tao->trust, bnk->min_radius); 462 tao->trust = PetscMin(tao->trust, bnk->max_radius); 463 } else { 464 /* The direction was bad; set radius to default value and re-solve 465 the trust-region subproblem to get a direction */ 466 tao->trust = tao->trust0; 467 468 /* Modify the radius if it is too large or small */ 469 tao->trust = PetscMax(tao->trust, bnk->min_radius); 470 tao->trust = PetscMin(tao->trust, bnk->max_radius); 471 472 PetscCall(KSPCGSetRadius(tao->ksp, tao->trust)); 473 PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive)); 474 PetscCall(KSPGetIterationNumber(tao->ksp, &kspits)); 475 tao->ksp_its += kspits; 476 tao->ksp_tot_its += kspits; 477 PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm)); 478 479 PetscCheck(bnk->dnorm != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero"); 480 } 481 } 482 } 483 /* Restore sub vectors back */ 484 if (bnk->active_idx) { 485 PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 486 PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 487 } 488 /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 489 PetscCall(VecScale(tao->stepdirection, -1.0)); 490 PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 491 492 /* Record convergence reasons */ 493 PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason)); 494 if (KSP_CONVERGED_ATOL == *ksp_reason) { 495 ++bnk->ksp_atol; 496 } else if (KSP_CONVERGED_RTOL == *ksp_reason) { 497 ++bnk->ksp_rtol; 498 } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) { 499 ++bnk->ksp_ctol; 500 } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) { 501 ++bnk->ksp_negc; 502 } else if (KSP_DIVERGED_DTOL == *ksp_reason) { 503 ++bnk->ksp_dtol; 504 } else if (KSP_DIVERGED_ITS == *ksp_reason) { 505 ++bnk->ksp_iter; 506 } else { 507 ++bnk->ksp_othr; 508 } 509 510 /* Make sure the BFGS preconditioner is healthy */ 511 if (bnk->M) { 512 PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates)); 513 if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) { 514 /* Preconditioner is numerically indefinite; reset the approximation. */ 515 PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 516 PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 517 } 518 } 519 *step_type = BNK_NEWTON; 520 PetscFunctionReturn(0); 521 } 522 523 /*------------------------------------------------------------*/ 524 525 /* Routine for recomputing the predicted reduction for a given step vector */ 526 527 PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered) { 528 TAO_BNK *bnk = (TAO_BNK *)tao->data; 529 530 PetscFunctionBegin; 531 /* Extract subvectors associated with the inactive set */ 532 if (bnk->active_idx) { 533 PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 534 PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 535 PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 536 } else { 537 bnk->X_inactive = tao->stepdirection; 538 bnk->inactive_work = bnk->Xwork; 539 bnk->G_inactive = bnk->Gwork; 540 } 541 /* Recompute the predicted decrease based on the quadratic model */ 542 PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work)); 543 PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive)); 544 PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered)); 545 /* Restore the sub vectors */ 546 if (bnk->active_idx) { 547 PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive)); 548 PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work)); 549 PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive)); 550 } 551 PetscFunctionReturn(0); 552 } 553 554 /*------------------------------------------------------------*/ 555 556 /* Routine for ensuring that the Newton step is a descent direction. 557 558 The step direction falls back onto BFGS, scaled gradient and gradient steps 559 in the event that the Newton step fails the test. 560 */ 561 562 PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType) { 563 TAO_BNK *bnk = (TAO_BNK *)tao->data; 564 PetscReal gdx, e_min; 565 PetscInt bfgsUpdates; 566 567 PetscFunctionBegin; 568 switch (*stepType) { 569 case BNK_NEWTON: 570 PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx)); 571 if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 572 /* Newton step is not descent or direction produced Inf or NaN 573 Update the perturbation for next time */ 574 if (bnk->pert <= 0.0) { 575 PetscBool is_gltr; 576 577 /* Initialize the perturbation */ 578 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 579 PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr)); 580 if (is_gltr) { 581 PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min)); 582 bnk->pert = PetscMax(bnk->pert, -e_min); 583 } 584 } else { 585 /* Increase the perturbation */ 586 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 587 } 588 589 if (!bnk->M) { 590 /* We don't have the bfgs matrix around and updated 591 Must use gradient direction in this case */ 592 PetscCall(VecCopy(tao->gradient, tao->stepdirection)); 593 *stepType = BNK_GRADIENT; 594 } else { 595 /* Attempt to use the BFGS direction */ 596 PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 597 598 /* Check for success (descent direction) 599 NOTE: Negative gdx here means not a descent direction because 600 the fall-back step is missing a negative sign. */ 601 PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 602 if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 603 /* BFGS direction is not descent or direction produced not a number 604 We can assert bfgsUpdates > 1 in this case because 605 the first solve produces the scaled gradient direction, 606 which is guaranteed to be descent */ 607 608 /* Use steepest descent direction (scaled) */ 609 PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 610 PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 611 PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 612 613 *stepType = BNK_SCALED_GRADIENT; 614 } else { 615 PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates)); 616 if (1 == bfgsUpdates) { 617 /* The first BFGS direction is always the scaled gradient */ 618 *stepType = BNK_SCALED_GRADIENT; 619 } else { 620 *stepType = BNK_BFGS; 621 } 622 } 623 } 624 /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 625 PetscCall(VecScale(tao->stepdirection, -1.0)); 626 PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 627 } else { 628 /* Computed Newton step is descent */ 629 switch (ksp_reason) { 630 case KSP_DIVERGED_NANORINF: 631 case KSP_DIVERGED_BREAKDOWN: 632 case KSP_DIVERGED_INDEFINITE_MAT: 633 case KSP_DIVERGED_INDEFINITE_PC: 634 case KSP_CONVERGED_CG_NEG_CURVE: 635 /* Matrix or preconditioner is indefinite; increase perturbation */ 636 if (bnk->pert <= 0.0) { 637 PetscBool is_gltr; 638 639 /* Initialize the perturbation */ 640 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 641 PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr)); 642 if (is_gltr) { 643 PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min)); 644 bnk->pert = PetscMax(bnk->pert, -e_min); 645 } 646 } else { 647 /* Increase the perturbation */ 648 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 649 } 650 break; 651 652 default: 653 /* Newton step computation is good; decrease perturbation */ 654 bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 655 if (bnk->pert < bnk->pmin) { bnk->pert = 0.0; } 656 break; 657 } 658 *stepType = BNK_NEWTON; 659 } 660 break; 661 662 case BNK_BFGS: 663 /* Check for success (descent direction) */ 664 PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx)); 665 if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) { 666 /* Step is not descent or solve was not successful 667 Use steepest descent direction (scaled) */ 668 PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 669 PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 670 PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection)); 671 PetscCall(VecScale(tao->stepdirection, -1.0)); 672 PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 673 *stepType = BNK_SCALED_GRADIENT; 674 } else { 675 *stepType = BNK_BFGS; 676 } 677 break; 678 679 case BNK_SCALED_GRADIENT: break; 680 681 default: break; 682 } 683 684 PetscFunctionReturn(0); 685 } 686 687 /*------------------------------------------------------------*/ 688 689 /* Routine for performing a bound-projected More-Thuente line search. 690 691 Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 692 Newton step does not produce a valid step length. 693 */ 694 695 PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) { 696 TAO_BNK *bnk = (TAO_BNK *)tao->data; 697 TaoLineSearchConvergedReason ls_reason; 698 PetscReal e_min, gdx; 699 PetscInt bfgsUpdates; 700 701 PetscFunctionBegin; 702 /* Perform the linesearch */ 703 PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason)); 704 PetscCall(TaoAddLineSearchCounts(tao)); 705 706 while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) { 707 /* Linesearch failed, revert solution */ 708 bnk->f = bnk->fold; 709 PetscCall(VecCopy(bnk->Xold, tao->solution)); 710 PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient)); 711 712 switch (*stepType) { 713 case BNK_NEWTON: 714 /* Failed to obtain acceptable iterate with Newton step 715 Update the perturbation for next time */ 716 if (bnk->pert <= 0.0) { 717 PetscBool is_gltr; 718 719 /* Initialize the perturbation */ 720 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 721 PetscCall(PetscObjectTypeCompare((PetscObject)(tao->ksp), KSPGLTR, &is_gltr)); 722 if (is_gltr) { 723 PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min)); 724 bnk->pert = PetscMax(bnk->pert, -e_min); 725 } 726 } else { 727 /* Increase the perturbation */ 728 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 729 } 730 731 if (!bnk->M) { 732 /* We don't have the bfgs matrix around and being updated 733 Must use gradient direction in this case */ 734 PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection)); 735 *stepType = BNK_GRADIENT; 736 } else { 737 /* Attempt to use the BFGS direction */ 738 PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 739 /* Check for success (descent direction) 740 NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 741 PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 742 if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) { 743 /* BFGS direction is not descent or direction produced not a number 744 We can assert bfgsUpdates > 1 in this case 745 Use steepest descent direction (scaled) */ 746 PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 747 PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 748 PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 749 750 bfgsUpdates = 1; 751 *stepType = BNK_SCALED_GRADIENT; 752 } else { 753 PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates)); 754 if (1 == bfgsUpdates) { 755 /* The first BFGS direction is always the scaled gradient */ 756 *stepType = BNK_SCALED_GRADIENT; 757 } else { 758 *stepType = BNK_BFGS; 759 } 760 } 761 } 762 break; 763 764 case BNK_BFGS: 765 /* Can only enter if pc_type == BNK_PC_BFGS 766 Failed to obtain acceptable iterate with BFGS step 767 Attempt to use the scaled gradient direction */ 768 PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE)); 769 PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient)); 770 PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection)); 771 772 bfgsUpdates = 1; 773 *stepType = BNK_SCALED_GRADIENT; 774 break; 775 } 776 /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 777 PetscCall(VecScale(tao->stepdirection, -1.0)); 778 PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection)); 779 780 /* Perform one last line search with the fall-back step */ 781 PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason)); 782 PetscCall(TaoAddLineSearchCounts(tao)); 783 } 784 *reason = ls_reason; 785 PetscFunctionReturn(0); 786 } 787 788 /*------------------------------------------------------------*/ 789 790 /* Routine for updating the trust radius. 791 792 Function features three different update methods: 793 1) Line-search step length based 794 2) Predicted decrease on the CG quadratic model 795 3) Interpolation 796 */ 797 798 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept) { 799 TAO_BNK *bnk = (TAO_BNK *)tao->data; 800 801 PetscReal step, kappa; 802 PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 803 804 PetscFunctionBegin; 805 /* Update trust region radius */ 806 *accept = PETSC_FALSE; 807 switch (updateType) { 808 case BNK_UPDATE_STEP: 809 *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 810 if (stepType == BNK_NEWTON) { 811 PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step)); 812 if (step < bnk->nu1) { 813 /* Very bad step taken; reduce radius */ 814 tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 815 } else if (step < bnk->nu2) { 816 /* Reasonably bad step taken; reduce radius */ 817 tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 818 } else if (step < bnk->nu3) { 819 /* Reasonable step was taken; leave radius alone */ 820 if (bnk->omega3 < 1.0) { 821 tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 822 } else if (bnk->omega3 > 1.0) { 823 tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 824 } 825 } else if (step < bnk->nu4) { 826 /* Full step taken; increase the radius */ 827 tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 828 } else { 829 /* More than full step taken; increase the radius */ 830 tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 831 } 832 } else { 833 /* Newton step was not good; reduce the radius */ 834 tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 835 } 836 break; 837 838 case BNK_UPDATE_REDUCTION: 839 if (stepType == BNK_NEWTON) { 840 if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) { 841 /* The predicted reduction has the wrong sign. This cannot 842 happen in infinite precision arithmetic. Step should 843 be rejected! */ 844 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 845 } else { 846 if (PetscIsInfOrNanReal(actred)) { 847 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 848 } else { 849 if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon)) { 850 kappa = 1.0; 851 } else { 852 kappa = actred / prered; 853 } 854 /* Accept or reject the step and update radius */ 855 if (kappa < bnk->eta1) { 856 /* Reject the step */ 857 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 858 } else { 859 /* Accept the step */ 860 *accept = PETSC_TRUE; 861 /* Update the trust region radius only if the computed step is at the trust radius boundary */ 862 if (bnk->dnorm == tao->trust) { 863 if (kappa < bnk->eta2) { 864 /* Marginal bad step */ 865 tao->trust = bnk->alpha2 * tao->trust; 866 } else if (kappa < bnk->eta3) { 867 /* Reasonable step */ 868 tao->trust = bnk->alpha3 * tao->trust; 869 } else if (kappa < bnk->eta4) { 870 /* Good step */ 871 tao->trust = bnk->alpha4 * tao->trust; 872 } else { 873 /* Very good step */ 874 tao->trust = bnk->alpha5 * tao->trust; 875 } 876 } 877 } 878 } 879 } 880 } else { 881 /* Newton step was not good; reduce the radius */ 882 tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 883 } 884 break; 885 886 default: 887 if (stepType == BNK_NEWTON) { 888 if (prered < 0.0) { 889 /* The predicted reduction has the wrong sign. This cannot */ 890 /* happen in infinite precision arithmetic. Step should */ 891 /* be rejected! */ 892 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 893 } else { 894 if (PetscIsInfOrNanReal(actred)) { 895 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 896 } else { 897 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 898 kappa = 1.0; 899 } else { 900 kappa = actred / prered; 901 } 902 903 PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 904 tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 905 tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 906 tau_min = PetscMin(tau_1, tau_2); 907 tau_max = PetscMax(tau_1, tau_2); 908 909 if (kappa >= 1.0 - bnk->mu1) { 910 /* Great agreement */ 911 *accept = PETSC_TRUE; 912 if (tau_max < 1.0) { 913 tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 914 } else if (tau_max > bnk->gamma4) { 915 tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 916 } else { 917 tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 918 } 919 } else if (kappa >= 1.0 - bnk->mu2) { 920 /* Good agreement */ 921 *accept = PETSC_TRUE; 922 if (tau_max < bnk->gamma2) { 923 tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 924 } else if (tau_max > bnk->gamma3) { 925 tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 926 } else if (tau_max < 1.0) { 927 tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 928 } else { 929 tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 930 } 931 } else { 932 /* Not good agreement */ 933 if (tau_min > 1.0) { 934 tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 935 } else if (tau_max < bnk->gamma1) { 936 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 937 } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 938 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 939 } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 940 tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 941 } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 942 tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 943 } else { 944 tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 945 } 946 } 947 } 948 } 949 } else { 950 /* Newton step was not good; reduce the radius */ 951 tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 952 } 953 break; 954 } 955 /* Make sure the radius does not violate min and max settings */ 956 tao->trust = PetscMin(tao->trust, bnk->max_radius); 957 tao->trust = PetscMax(tao->trust, bnk->min_radius); 958 PetscFunctionReturn(0); 959 } 960 961 /* ---------------------------------------------------------- */ 962 963 PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType) { 964 TAO_BNK *bnk = (TAO_BNK *)tao->data; 965 966 PetscFunctionBegin; 967 switch (stepType) { 968 case BNK_NEWTON: ++bnk->newt; break; 969 case BNK_BFGS: ++bnk->bfgs; break; 970 case BNK_SCALED_GRADIENT: ++bnk->sgrad; break; 971 case BNK_GRADIENT: ++bnk->grad; break; 972 default: break; 973 } 974 PetscFunctionReturn(0); 975 } 976 977 /* ---------------------------------------------------------- */ 978 979 PetscErrorCode TaoSetUp_BNK(Tao tao) { 980 TAO_BNK *bnk = (TAO_BNK *)tao->data; 981 PetscInt i; 982 983 PetscFunctionBegin; 984 if (!tao->gradient) { PetscCall(VecDuplicate(tao->solution, &tao->gradient)); } 985 if (!tao->stepdirection) { PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); } 986 if (!bnk->W) { PetscCall(VecDuplicate(tao->solution, &bnk->W)); } 987 if (!bnk->Xold) { PetscCall(VecDuplicate(tao->solution, &bnk->Xold)); } 988 if (!bnk->Gold) { PetscCall(VecDuplicate(tao->solution, &bnk->Gold)); } 989 if (!bnk->Xwork) { PetscCall(VecDuplicate(tao->solution, &bnk->Xwork)); } 990 if (!bnk->Gwork) { PetscCall(VecDuplicate(tao->solution, &bnk->Gwork)); } 991 if (!bnk->unprojected_gradient) { PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient)); } 992 if (!bnk->unprojected_gradient_old) { PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient_old)); } 993 if (!bnk->Diag_min) { PetscCall(VecDuplicate(tao->solution, &bnk->Diag_min)); } 994 if (!bnk->Diag_max) { PetscCall(VecDuplicate(tao->solution, &bnk->Diag_max)); } 995 if (bnk->max_cg_its > 0) { 996 /* Ensure that the important common vectors are shared between BNK and embedded BNCG */ 997 bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data; 998 PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old))); 999 PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old)); 1000 bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old; 1001 PetscCall(PetscObjectReference((PetscObject)(bnk->unprojected_gradient))); 1002 PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient)); 1003 bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient; 1004 PetscCall(PetscObjectReference((PetscObject)(bnk->Gold))); 1005 PetscCall(VecDestroy(&bnk->bncg_ctx->G_old)); 1006 bnk->bncg_ctx->G_old = bnk->Gold; 1007 PetscCall(PetscObjectReference((PetscObject)(tao->gradient))); 1008 PetscCall(VecDestroy(&bnk->bncg->gradient)); 1009 bnk->bncg->gradient = tao->gradient; 1010 PetscCall(PetscObjectReference((PetscObject)(tao->stepdirection))); 1011 PetscCall(VecDestroy(&bnk->bncg->stepdirection)); 1012 bnk->bncg->stepdirection = tao->stepdirection; 1013 PetscCall(TaoSetSolution(bnk->bncg, tao->solution)); 1014 /* Copy over some settings from BNK into BNCG */ 1015 PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its)); 1016 PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol)); 1017 PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin)); 1018 PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP)); 1019 PetscCall(TaoSetObjective(bnk->bncg, tao->ops->computeobjective, tao->user_objP)); 1020 PetscCall(TaoSetGradient(bnk->bncg, NULL, tao->ops->computegradient, tao->user_gradP)); 1021 PetscCall(TaoSetObjectiveAndGradient(bnk->bncg, NULL, tao->ops->computeobjectiveandgradient, tao->user_objgradP)); 1022 PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg))); 1023 for (i = 0; i < tao->numbermonitors; ++i) { 1024 PetscCall(TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i])); 1025 PetscCall(PetscObjectReference((PetscObject)(tao->monitorcontext[i]))); 1026 } 1027 } 1028 bnk->X_inactive = NULL; 1029 bnk->G_inactive = NULL; 1030 bnk->inactive_work = NULL; 1031 bnk->active_work = NULL; 1032 bnk->inactive_idx = NULL; 1033 bnk->active_idx = NULL; 1034 bnk->active_lower = NULL; 1035 bnk->active_upper = NULL; 1036 bnk->active_fixed = NULL; 1037 bnk->M = NULL; 1038 bnk->H_inactive = NULL; 1039 bnk->Hpre_inactive = NULL; 1040 PetscFunctionReturn(0); 1041 } 1042 1043 /*------------------------------------------------------------*/ 1044 1045 PetscErrorCode TaoDestroy_BNK(Tao tao) { 1046 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1047 1048 PetscFunctionBegin; 1049 PetscCall(VecDestroy(&bnk->W)); 1050 PetscCall(VecDestroy(&bnk->Xold)); 1051 PetscCall(VecDestroy(&bnk->Gold)); 1052 PetscCall(VecDestroy(&bnk->Xwork)); 1053 PetscCall(VecDestroy(&bnk->Gwork)); 1054 PetscCall(VecDestroy(&bnk->unprojected_gradient)); 1055 PetscCall(VecDestroy(&bnk->unprojected_gradient_old)); 1056 PetscCall(VecDestroy(&bnk->Diag_min)); 1057 PetscCall(VecDestroy(&bnk->Diag_max)); 1058 PetscCall(ISDestroy(&bnk->active_lower)); 1059 PetscCall(ISDestroy(&bnk->active_upper)); 1060 PetscCall(ISDestroy(&bnk->active_fixed)); 1061 PetscCall(ISDestroy(&bnk->active_idx)); 1062 PetscCall(ISDestroy(&bnk->inactive_idx)); 1063 PetscCall(MatDestroy(&bnk->Hpre_inactive)); 1064 PetscCall(MatDestroy(&bnk->H_inactive)); 1065 PetscCall(TaoDestroy(&bnk->bncg)); 1066 PetscCall(KSPDestroy(&tao->ksp)); 1067 PetscCall(PetscFree(tao->data)); 1068 PetscFunctionReturn(0); 1069 } 1070 1071 /*------------------------------------------------------------*/ 1072 1073 PetscErrorCode TaoSetFromOptions_BNK(Tao tao, PetscOptionItems *PetscOptionsObject) { 1074 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1075 1076 PetscFunctionBegin; 1077 PetscOptionsHeadBegin(PetscOptionsObject, "Newton-Krylov method for bound constrained optimization"); 1078 PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL)); 1079 PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL)); 1080 PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL)); 1081 PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval, NULL)); 1082 PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin, NULL)); 1083 PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax, NULL)); 1084 PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac, NULL)); 1085 PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin, NULL)); 1086 PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax, NULL)); 1087 PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac, NULL)); 1088 PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac, NULL)); 1089 PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac, NULL)); 1090 PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac, NULL)); 1091 PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1, NULL)); 1092 PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2, NULL)); 1093 PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3, NULL)); 1094 PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4, NULL)); 1095 PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1, NULL)); 1096 PetscCall(PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2, NULL)); 1097 PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3, NULL)); 1098 PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4, NULL)); 1099 PetscCall(PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5, NULL)); 1100 PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1, NULL)); 1101 PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2, NULL)); 1102 PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3, NULL)); 1103 PetscCall(PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4, NULL)); 1104 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)); 1105 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)); 1106 PetscCall(PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3, NULL)); 1107 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)); 1108 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)); 1109 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)); 1110 PetscCall(PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i, NULL)); 1111 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)); 1112 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)); 1113 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)); 1114 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)); 1115 PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i, NULL)); 1116 PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1, NULL)); 1117 PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2, NULL)); 1118 PetscCall(PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1, NULL)); 1119 PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2, NULL)); 1120 PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3, NULL)); 1121 PetscCall(PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4, NULL)); 1122 PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta, NULL)); 1123 PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius, NULL)); 1124 PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius, NULL)); 1125 PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL)); 1126 PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL)); 1127 PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL)); 1128 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)); 1129 PetscOptionsHeadEnd(); 1130 1131 PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)(tao))->prefix)); 1132 PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_cg_")); 1133 PetscCall(TaoSetFromOptions(bnk->bncg)); 1134 1135 PetscCall(KSPSetOptionsPrefix(tao->ksp, ((PetscObject)(tao))->prefix)); 1136 PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_bnk_")); 1137 PetscCall(KSPSetFromOptions(tao->ksp)); 1138 PetscFunctionReturn(0); 1139 } 1140 1141 /*------------------------------------------------------------*/ 1142 1143 PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) { 1144 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1145 PetscInt nrejects; 1146 PetscBool isascii; 1147 1148 PetscFunctionBegin; 1149 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 1150 if (isascii) { 1151 PetscCall(PetscViewerASCIIPushTab(viewer)); 1152 if (bnk->M) { 1153 PetscCall(MatLMVMGetRejectCount(bnk->M, &nrejects)); 1154 PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n", nrejects)); 1155 } 1156 PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its)); 1157 PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt)); 1158 if (bnk->M) { PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs)); } 1159 PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad)); 1160 PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad)); 1161 PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n")); 1162 PetscCall(PetscViewerASCIIPrintf(viewer, " atol: %" PetscInt_FMT "\n", bnk->ksp_atol)); 1163 PetscCall(PetscViewerASCIIPrintf(viewer, " rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol)); 1164 PetscCall(PetscViewerASCIIPrintf(viewer, " ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol)); 1165 PetscCall(PetscViewerASCIIPrintf(viewer, " negc: %" PetscInt_FMT "\n", bnk->ksp_negc)); 1166 PetscCall(PetscViewerASCIIPrintf(viewer, " dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol)); 1167 PetscCall(PetscViewerASCIIPrintf(viewer, " iter: %" PetscInt_FMT "\n", bnk->ksp_iter)); 1168 PetscCall(PetscViewerASCIIPrintf(viewer, " othr: %" PetscInt_FMT "\n", bnk->ksp_othr)); 1169 PetscCall(PetscViewerASCIIPopTab(viewer)); 1170 } 1171 PetscFunctionReturn(0); 1172 } 1173 1174 /* ---------------------------------------------------------- */ 1175 1176 /*MC 1177 TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 1178 At each iteration, the BNK methods solve the symmetric 1179 system of equations to obtain the step diretion dk: 1180 Hk dk = -gk 1181 for free variables only. The step can be globalized either through 1182 trust-region methods, or a line search, or a heuristic mixture of both. 1183 1184 Options Database Keys: 1185 + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop 1186 . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation") 1187 . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation") 1188 . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas") 1189 . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas) 1190 . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas) 1191 . -tao_bnk_sval - (developer) Hessian perturbation starting value 1192 . -tao_bnk_imin - (developer) minimum initial Hessian perturbation 1193 . -tao_bnk_imax - (developer) maximum initial Hessian perturbation 1194 . -tao_bnk_pmin - (developer) minimum Hessian perturbation 1195 . -tao_bnk_pmax - (developer) aximum Hessian perturbation 1196 . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor 1197 . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor 1198 . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation 1199 . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation 1200 . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation 1201 . -tao_bnk_eta1 - (developer) threshold for rejecting step (-update_type reduction) 1202 . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-update_type reduction) 1203 . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-update_type reduction) 1204 . -tao_bnk_eta4 - (developer) threshold for accepting good step (-update_type reduction) 1205 . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction) 1206 . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction) 1207 . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction) 1208 . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction) 1209 . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction) 1210 . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction) 1211 . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-update_type interpolation) 1212 . -tao_bnk_mu2 - (developer) threshold for accepting good step (-update_type interpolation) 1213 . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation) 1214 . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation) 1215 . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation) 1216 . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation) 1217 . -tao_bnk_theta - (developer) trust region interpolation factor (-update_type interpolation) 1218 . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-update_type step) 1219 . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-update_type step) 1220 . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-update_type step) 1221 . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-update_type step) 1222 . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step) 1223 . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-update_type step) 1224 . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-update_type step) 1225 . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-update_type step) 1226 . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-update_type step) 1227 . -tao_bnk_mu1_i - (developer) threshold for accepting very good step (-init_type interpolation) 1228 . -tao_bnk_mu2_i - (developer) threshold for accepting good step (-init_type interpolation) 1229 . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation) 1230 . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation) 1231 . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation) 1232 . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation) 1233 - -tao_bnk_theta_i - (developer) trust region interpolation factor (-init_type interpolation) 1234 1235 Level: beginner 1236 M*/ 1237 1238 PetscErrorCode TaoCreate_BNK(Tao tao) { 1239 TAO_BNK *bnk; 1240 PC pc; 1241 1242 PetscFunctionBegin; 1243 PetscCall(PetscNewLog(tao, &bnk)); 1244 1245 tao->ops->setup = TaoSetUp_BNK; 1246 tao->ops->view = TaoView_BNK; 1247 tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1248 tao->ops->destroy = TaoDestroy_BNK; 1249 1250 /* Override default settings (unless already changed) */ 1251 if (!tao->max_it_changed) tao->max_it = 50; 1252 if (!tao->trust0_changed) tao->trust0 = 100.0; 1253 1254 tao->data = (void *)bnk; 1255 1256 /* Hessian shifting parameters */ 1257 bnk->computehessian = TaoBNKComputeHessian; 1258 bnk->computestep = TaoBNKComputeStep; 1259 1260 bnk->sval = 0.0; 1261 bnk->imin = 1.0e-4; 1262 bnk->imax = 1.0e+2; 1263 bnk->imfac = 1.0e-1; 1264 1265 bnk->pmin = 1.0e-12; 1266 bnk->pmax = 1.0e+2; 1267 bnk->pgfac = 1.0e+1; 1268 bnk->psfac = 4.0e-1; 1269 bnk->pmgfac = 1.0e-1; 1270 bnk->pmsfac = 1.0e-1; 1271 1272 /* Default values for trust-region radius update based on steplength */ 1273 bnk->nu1 = 0.25; 1274 bnk->nu2 = 0.50; 1275 bnk->nu3 = 1.00; 1276 bnk->nu4 = 1.25; 1277 1278 bnk->omega1 = 0.25; 1279 bnk->omega2 = 0.50; 1280 bnk->omega3 = 1.00; 1281 bnk->omega4 = 2.00; 1282 bnk->omega5 = 4.00; 1283 1284 /* Default values for trust-region radius update based on reduction */ 1285 bnk->eta1 = 1.0e-4; 1286 bnk->eta2 = 0.25; 1287 bnk->eta3 = 0.50; 1288 bnk->eta4 = 0.90; 1289 1290 bnk->alpha1 = 0.25; 1291 bnk->alpha2 = 0.50; 1292 bnk->alpha3 = 1.00; 1293 bnk->alpha4 = 2.00; 1294 bnk->alpha5 = 4.00; 1295 1296 /* Default values for trust-region radius update based on interpolation */ 1297 bnk->mu1 = 0.10; 1298 bnk->mu2 = 0.50; 1299 1300 bnk->gamma1 = 0.25; 1301 bnk->gamma2 = 0.50; 1302 bnk->gamma3 = 2.00; 1303 bnk->gamma4 = 4.00; 1304 1305 bnk->theta = 0.05; 1306 1307 /* Default values for trust region initialization based on interpolation */ 1308 bnk->mu1_i = 0.35; 1309 bnk->mu2_i = 0.50; 1310 1311 bnk->gamma1_i = 0.0625; 1312 bnk->gamma2_i = 0.5; 1313 bnk->gamma3_i = 2.0; 1314 bnk->gamma4_i = 5.0; 1315 1316 bnk->theta_i = 0.25; 1317 1318 /* Remaining parameters */ 1319 bnk->max_cg_its = 0; 1320 bnk->min_radius = 1.0e-10; 1321 bnk->max_radius = 1.0e10; 1322 bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0); 1323 bnk->as_tol = 1.0e-3; 1324 bnk->as_step = 1.0e-3; 1325 bnk->dmin = 1.0e-6; 1326 bnk->dmax = 1.0e6; 1327 1328 bnk->M = NULL; 1329 bnk->bfgs_pre = NULL; 1330 bnk->init_type = BNK_INIT_INTERPOLATION; 1331 bnk->update_type = BNK_UPDATE_REDUCTION; 1332 bnk->as_type = BNK_AS_BERTSEKAS; 1333 1334 /* Create the embedded BNCG solver */ 1335 PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg)); 1336 PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1)); 1337 PetscCall(TaoSetType(bnk->bncg, TAOBNCG)); 1338 1339 /* Create the line search */ 1340 PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 1341 PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 1342 PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT)); 1343 PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao)); 1344 1345 /* Set linear solver to default for symmetric matrices */ 1346 PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 1347 PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 1348 PetscCall(KSPSetType(tao->ksp, KSPSTCG)); 1349 PetscCall(KSPGetPC(tao->ksp, &pc)); 1350 PetscCall(PCSetType(pc, PCLMVM)); 1351 PetscFunctionReturn(0); 1352 } 1353