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