1 #include <petsctaolinesearch.h> 2 #include <../src/tao/bound/impls/bnk/bnk.h> 3 4 #include <petscksp.h> 5 6 /* Routine for BFGS preconditioner */ 7 8 PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 9 { 10 PetscErrorCode ierr; 11 Mat M; 12 13 PetscFunctionBegin; 14 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15 PetscValidHeaderSpecific(b,VEC_CLASSID,2); 16 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 17 ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 18 ierr = MatLMVMSolveInactive(M, b, x);CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 22 /*------------------------------------------------------------*/ 23 24 /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */ 25 26 PetscErrorCode TaoBNKInitialize(Tao tao) 27 { 28 PetscErrorCode ierr; 29 TAO_BNK *bnk = (TAO_BNK *)tao->data; 30 KSPType ksp_type; 31 PC pc; 32 33 PetscReal fmin, ftrial, prered, actred, kappa, sigma, dnorm; 34 PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 35 PetscReal delta, step = 1.0; 36 37 PetscInt n,N,needH = 1; 38 39 PetscInt i_max = 5; 40 PetscInt j_max = 1; 41 PetscInt i, j; 42 43 PetscFunctionBegin; 44 /* Number of times ksp stopped because of these reasons */ 45 bnk->ksp_atol = 0; 46 bnk->ksp_rtol = 0; 47 bnk->ksp_dtol = 0; 48 bnk->ksp_ctol = 0; 49 bnk->ksp_negc = 0; 50 bnk->ksp_iter = 0; 51 bnk->ksp_othr = 0; 52 53 /* Initialize the Hessian perturbation */ 54 bnk->pert = bnk->sval; 55 56 /* Initialize trust-region radius when using nash, stcg, or gltr 57 Command automatically ignored for other methods 58 Will be reset during the first iteration 59 */ 60 ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 61 ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 62 ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 63 ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 64 65 ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 66 67 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 68 if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative"); 69 tao->trust = tao->trust0; 70 tao->trust = PetscMax(tao->trust, bnk->min_radius); 71 tao->trust = PetscMin(tao->trust, bnk->max_radius); 72 } 73 74 /* Get vectors we will need */ 75 if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) { 76 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 77 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 78 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr); 79 ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr); 80 } 81 82 /* create vectors for the limited memory preconditioner */ 83 if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) { 84 if (!bnk->Diag) { 85 ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr); 86 } 87 } 88 89 /* Modify the preconditioner to use the bfgs approximation */ 90 ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 91 switch(bnk->pc_type) { 92 case BNK_PC_NONE: 93 ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 94 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 95 break; 96 97 case BNK_PC_AHESS: 98 ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 99 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 100 ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 101 break; 102 103 case BNK_PC_BFGS: 104 ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 105 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 106 ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 107 ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr); 108 ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 109 break; 110 111 default: 112 /* Use the pc method set by pc_type */ 113 break; 114 } 115 116 /* Initialize trust-region radius. The initialization is only performed 117 when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 118 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 119 switch(bnk->init_type) { 120 case BNK_INIT_CONSTANT: 121 /* Use the initial radius specified */ 122 break; 123 124 case BNK_INIT_INTERPOLATION: 125 /* Use the initial radius specified */ 126 max_radius = 0.0; 127 128 for (j = 0; j < j_max; ++j) { 129 fmin = bnk->f; 130 sigma = 0.0; 131 132 if (needH) { 133 ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 134 needH = 0; 135 } 136 137 for (i = 0; i < i_max; ++i) { 138 ierr = VecCopy(tao->solution,bnk->W);CHKERRQ(ierr); 139 ierr = VecAXPY(bnk->W,-tao->trust/bnk->gnorm,tao->gradient);CHKERRQ(ierr); 140 ierr = VecMedian(tao->XL, bnk->W, tao->XU, bnk->W);CHKERRQ(ierr); 141 ierr = VecNorm(bnk->W, NORM_2, &dnorm);CHKERRQ(ierr); 142 if (dnorm != tao->trust/bnk->gnorm) tao->trust = dnorm; 143 ierr = TaoComputeObjective(tao, bnk->W, &ftrial);CHKERRQ(ierr); 144 if (PetscIsInfOrNanReal(ftrial)) { 145 tau = bnk->gamma1_i; 146 } else { 147 if (ftrial < fmin) { 148 fmin = ftrial; 149 sigma = -tao->trust / bnk->gnorm; 150 } 151 152 ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 153 ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 154 155 prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 156 actred = bnk->f - ftrial; 157 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 158 kappa = 1.0; 159 } else { 160 kappa = actred / prered; 161 } 162 163 tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 164 tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 165 tau_min = PetscMin(tau_1, tau_2); 166 tau_max = PetscMax(tau_1, tau_2); 167 168 if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) { 169 /* Great agreement */ 170 max_radius = PetscMax(max_radius, tao->trust); 171 172 if (tau_max < 1.0) { 173 tau = bnk->gamma3_i; 174 } else if (tau_max > bnk->gamma4_i) { 175 tau = bnk->gamma4_i; 176 } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) { 177 tau = tau_1; 178 } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) { 179 tau = tau_2; 180 } else { 181 tau = tau_max; 182 } 183 } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) { 184 /* Good agreement */ 185 max_radius = PetscMax(max_radius, tao->trust); 186 187 if (tau_max < bnk->gamma2_i) { 188 tau = bnk->gamma2_i; 189 } else if (tau_max > bnk->gamma3_i) { 190 tau = bnk->gamma3_i; 191 } else { 192 tau = tau_max; 193 } 194 } else { 195 /* Not good agreement */ 196 if (tau_min > 1.0) { 197 tau = bnk->gamma2_i; 198 } else if (tau_max < bnk->gamma1_i) { 199 tau = bnk->gamma1_i; 200 } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 201 tau = bnk->gamma1_i; 202 } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 203 tau = tau_1; 204 } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 205 tau = tau_2; 206 } else { 207 tau = tau_max; 208 } 209 } 210 } 211 tao->trust = tau * tao->trust; 212 } 213 214 if (fmin < bnk->f) { 215 bnk->f = fmin; 216 ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 217 ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 218 ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 219 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 220 221 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 222 if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 223 needH = 1; 224 225 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 226 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 227 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 228 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 229 } 230 } 231 tao->trust = PetscMax(tao->trust, max_radius); 232 233 /* Modify the radius if it is too large or small */ 234 tao->trust = PetscMax(tao->trust, bnk->min_radius); 235 tao->trust = PetscMin(tao->trust, bnk->max_radius); 236 break; 237 238 default: 239 /* Norm of the first direction will initialize radius */ 240 tao->trust = 0.0; 241 break; 242 } 243 } 244 245 /* Set initial scaling for the BFGS preconditioner 246 This step is done after computing the initial trust-region radius 247 since the function value may have decreased */ 248 if (BNK_PC_BFGS == bnk->pc_type) { 249 if (bnk->f != 0.0) { 250 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 251 } else { 252 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 253 } 254 ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 255 } 256 257 /* Set counter for gradient/reset steps*/ 258 bnk->newt = 0; 259 bnk->bfgs = 0; 260 bnk->sgrad = 0; 261 bnk->grad = 0; 262 PetscFunctionReturn(0); 263 } 264 265 /*------------------------------------------------------------*/ 266 267 /* Routine for computing the Newton step. 268 269 If the safeguard is enabled, the Newton step is verified to be a 270 descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled 271 gradient steps if/when necessary. 272 273 The function reports back on which type of step has ultimately been stored 274 under tao->stepdirection. 275 */ 276 277 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool safeguard, PetscInt *stepType) 278 { 279 PetscErrorCode ierr; 280 TAO_BNK *bnk = (TAO_BNK *)tao->data; 281 KSPConvergedReason ksp_reason; 282 283 PetscReal gdx, delta, e_min; 284 285 PetscInt bfgsUpdates = 0; 286 PetscInt kspits; 287 288 PetscFunctionBegin; 289 /* Shift the Hessian matrix */ 290 if (bnk->pert > 0) { 291 ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr); 292 if (tao->hessian != tao->hessian_pre) { 293 ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr); 294 } 295 } 296 297 /* Determine the inactive set */ 298 ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 299 ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr); 300 301 /* Prepare masked matrices for the inactive set */ 302 ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); 303 ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr); 304 ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 305 if (tao->hessian == tao->hessian_pre) { 306 bnk->Hpre_inactive = bnk->H_inactive; 307 } else { 308 ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr); 309 } 310 311 /* Solve the Newton system of equations */ 312 ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr); 313 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 314 ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 315 ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 316 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 317 tao->ksp_its+=kspits; 318 tao->ksp_tot_its+=kspits; 319 ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 320 321 if (0.0 == tao->trust) { 322 /* Radius was uninitialized; use the norm of the direction */ 323 if (bnk->dnorm > 0.0) { 324 tao->trust = bnk->dnorm; 325 326 /* Modify the radius if it is too large or small */ 327 tao->trust = PetscMax(tao->trust, bnk->min_radius); 328 tao->trust = PetscMin(tao->trust, bnk->max_radius); 329 } else { 330 /* The direction was bad; set radius to default value and re-solve 331 the trust-region subproblem to get a direction */ 332 tao->trust = tao->trust0; 333 334 /* Modify the radius if it is too large or small */ 335 tao->trust = PetscMax(tao->trust, bnk->min_radius); 336 tao->trust = PetscMin(tao->trust, bnk->max_radius); 337 338 ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 339 ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 340 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 341 tao->ksp_its+=kspits; 342 tao->ksp_tot_its+=kspits; 343 ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 344 345 if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 346 } 347 } 348 } else { 349 ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 350 ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 351 tao->ksp_its += kspits; 352 tao->ksp_tot_its+=kspits; 353 } 354 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 355 ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 356 357 /* Destroy masked matrices */ 358 if (bnk->H_inactive != bnk->Hpre_inactive) { 359 ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr); 360 } 361 ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 362 363 /* Make sure the BFGS preconditioner is healthy */ 364 if (bnk->pc_type == BNK_PC_BFGS) { 365 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 366 if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (bfgsUpdates > 1)) { 367 /* Preconditioner is numerically indefinite; reset the approximation. */ 368 if (bnk->f != 0.0) { 369 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 370 } else { 371 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 372 } 373 ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 374 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 375 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 376 bfgsUpdates = 1; 377 } 378 } 379 380 if (KSP_CONVERGED_ATOL == ksp_reason) { 381 ++bnk->ksp_atol; 382 } else if (KSP_CONVERGED_RTOL == ksp_reason) { 383 ++bnk->ksp_rtol; 384 } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 385 ++bnk->ksp_ctol; 386 } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 387 ++bnk->ksp_negc; 388 } else if (KSP_DIVERGED_DTOL == ksp_reason) { 389 ++bnk->ksp_dtol; 390 } else if (KSP_DIVERGED_ITS == ksp_reason) { 391 ++bnk->ksp_iter; 392 } else { 393 ++bnk->ksp_othr; 394 } 395 396 /* Check for success (descent direction) */ 397 if (safeguard) { 398 ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 399 if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 400 /* Newton step is not descent or direction produced Inf or NaN 401 Update the perturbation for next time */ 402 if (bnk->pert <= 0.0) { 403 /* Initialize the perturbation */ 404 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 405 if (bnk->is_gltr) { 406 ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 407 bnk->pert = PetscMax(bnk->pert, -e_min); 408 } 409 } else { 410 /* Increase the perturbation */ 411 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 412 } 413 414 if (BNK_PC_BFGS != bnk->pc_type) { 415 /* We don't have the bfgs matrix around and updated 416 Must use gradient direction in this case */ 417 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 418 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 419 ++bnk->grad; 420 *stepType = BNK_GRADIENT; 421 } else { 422 /* Attempt to use the BFGS direction */ 423 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 424 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 425 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 426 427 /* Check for success (descent direction) */ 428 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 429 if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 430 /* BFGS direction is not descent or direction produced not a number 431 We can assert bfgsUpdates > 1 in this case because 432 the first solve produces the scaled gradient direction, 433 which is guaranteed to be descent */ 434 435 /* Use steepest descent direction (scaled) */ 436 437 if (bnk->f != 0.0) { 438 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 439 } else { 440 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 441 } 442 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 443 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 444 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 445 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 446 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 447 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 448 449 bfgsUpdates = 1; 450 ++bnk->sgrad; 451 *stepType = BNK_SCALED_GRADIENT; 452 } else { 453 if (1 == bfgsUpdates) { 454 /* The first BFGS direction is always the scaled gradient */ 455 ++bnk->sgrad; 456 *stepType = BNK_SCALED_GRADIENT; 457 } else { 458 ++bnk->bfgs; 459 *stepType = BNK_BFGS; 460 } 461 } 462 } 463 } else { 464 /* Computed Newton step is descent */ 465 switch (ksp_reason) { 466 case KSP_DIVERGED_NANORINF: 467 case KSP_DIVERGED_BREAKDOWN: 468 case KSP_DIVERGED_INDEFINITE_MAT: 469 case KSP_DIVERGED_INDEFINITE_PC: 470 case KSP_CONVERGED_CG_NEG_CURVE: 471 /* Matrix or preconditioner is indefinite; increase perturbation */ 472 if (bnk->pert <= 0.0) { 473 /* Initialize the perturbation */ 474 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 475 if (bnk->is_gltr) { 476 ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 477 bnk->pert = PetscMax(bnk->pert, -e_min); 478 } 479 } else { 480 /* Increase the perturbation */ 481 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 482 } 483 break; 484 485 default: 486 /* Newton step computation is good; decrease perturbation */ 487 bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 488 if (bnk->pert < bnk->pmin) { 489 bnk->pert = 0.0; 490 } 491 break; 492 } 493 494 ++bnk->newt; 495 *stepType = BNK_NEWTON; 496 } 497 } else { 498 ++bnk->newt; 499 bnk->pert = 0.0; 500 *stepType = BNK_NEWTON; 501 } 502 PetscFunctionReturn(0); 503 } 504 505 /*------------------------------------------------------------*/ 506 507 /* Routine for performing a bound-projected More-Thuente line search. 508 509 Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 510 Newton step does not produce a valid step length. 511 */ 512 513 PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 514 { 515 TAO_BNK *bnk = (TAO_BNK *)tao->data; 516 PetscErrorCode ierr; 517 TaoLineSearchConvergedReason ls_reason; 518 519 PetscReal e_min, gdx, delta; 520 PetscInt bfgsUpdates; 521 522 PetscFunctionBegin; 523 /* Perform the linesearch */ 524 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 525 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 526 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 527 528 while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) { 529 /* Linesearch failed, revert solution */ 530 bnk->f = bnk->fold; 531 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 532 ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 533 ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 534 535 switch(stepType) { 536 case BNK_NEWTON: 537 /* Failed to obtain acceptable iterate with Newton 1step 538 Update the perturbation for next time */ 539 --bnk->newt; 540 if (bnk->pert <= 0.0) { 541 /* Initialize the perturbation */ 542 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 543 if (bnk->is_gltr) { 544 ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 545 bnk->pert = PetscMax(bnk->pert, -e_min); 546 } 547 } else { 548 /* Increase the perturbation */ 549 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 550 } 551 552 if (BNK_PC_BFGS != bnk->pc_type) { 553 /* We don't have the bfgs matrix around and being updated 554 Must use gradient direction in this case */ 555 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 556 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 557 ++bnk->grad; 558 stepType = BNK_GRADIENT; 559 } else { 560 /* Attempt to use the BFGS direction */ 561 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 562 ierr = VecBoundGradientProjection(tao->stepdirection, tao->solution, tao->XL, tao->XU, tao->stepdirection);CHKERRQ(ierr); 563 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 564 /* Check for success (descent direction) */ 565 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 566 if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 567 /* BFGS direction is not descent or direction produced not a number 568 We can assert bfgsUpdates > 1 in this case 569 Use steepest descent direction (scaled) */ 570 571 if (bnk->f != 0.0) { 572 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 573 } else { 574 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 575 } 576 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 577 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 578 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 579 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 580 ierr = VecBoundGradientProjection(tao->stepdirection, tao->solution, tao->XL, tao->XU, tao->stepdirection);CHKERRQ(ierr); 581 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 582 583 bfgsUpdates = 1; 584 ++bnk->sgrad; 585 stepType = BNK_SCALED_GRADIENT; 586 } else { 587 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 588 if (1 == bfgsUpdates) { 589 /* The first BFGS direction is always the scaled gradient */ 590 ++bnk->sgrad; 591 stepType = BNK_SCALED_GRADIENT; 592 } else { 593 ++bnk->bfgs; 594 stepType = BNK_BFGS; 595 } 596 } 597 } 598 break; 599 600 case BNK_BFGS: 601 /* Can only enter if pc_type == BNK_PC_BFGS 602 Failed to obtain acceptable iterate with BFGS step 603 Attempt to use the scaled gradient direction */ 604 605 if (bnk->f != 0.0) { 606 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 607 } else { 608 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 609 } 610 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 611 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 612 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 613 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 614 ierr = VecBoundGradientProjection(tao->stepdirection, tao->solution, tao->XL, tao->XU, tao->stepdirection);CHKERRQ(ierr); 615 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 616 617 bfgsUpdates = 1; 618 ++bnk->sgrad; 619 stepType = BNK_SCALED_GRADIENT; 620 break; 621 622 case BNK_SCALED_GRADIENT: 623 /* Can only enter if pc_type == BNK_PC_BFGS 624 The scaled gradient step did not produce a new iterate; 625 attemp to use the gradient direction. 626 Need to make sure we are not using a different diagonal scaling */ 627 628 ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 629 ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 630 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 631 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 632 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 633 ierr = VecBoundGradientProjection(tao->stepdirection, tao->solution, tao->XL, tao->XU, tao->stepdirection);CHKERRQ(ierr); 634 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 635 636 bfgsUpdates = 1; 637 ++bnk->grad; 638 stepType = BNK_GRADIENT; 639 break; 640 } 641 642 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 643 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 644 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 645 } 646 *reason = ls_reason; 647 PetscFunctionReturn(0); 648 } 649 650 /*------------------------------------------------------------*/ 651 652 /* Routine for updating the trust radius. 653 654 Function features three different update methods: 655 1) Line-search step length based 656 2) Predicted decrease on the CG quadratic model 657 3) Interpolation 658 */ 659 660 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt stepType, PetscBool *accept) 661 { 662 TAO_BNK *bnk = (TAO_BNK *)tao->data; 663 PetscErrorCode ierr; 664 665 PetscReal step, kappa; 666 PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 667 668 PetscFunctionBegin; 669 /* Update trust region radius */ 670 *accept = PETSC_FALSE; 671 switch(bnk->update_type) { 672 case BNK_UPDATE_STEP: 673 *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 674 if (stepType == BNK_NEWTON) { 675 ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr); 676 if (step < bnk->nu1) { 677 /* Very bad step taken; reduce radius */ 678 tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 679 } else if (step < bnk->nu2) { 680 /* Reasonably bad step taken; reduce radius */ 681 tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 682 } else if (step < bnk->nu3) { 683 /* Reasonable step was taken; leave radius alone */ 684 if (bnk->omega3 < 1.0) { 685 tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 686 } else if (bnk->omega3 > 1.0) { 687 tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 688 } 689 } else if (step < bnk->nu4) { 690 /* Full step taken; increase the radius */ 691 tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 692 } else { 693 /* More than full step taken; increase the radius */ 694 tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 695 } 696 } else { 697 /* Newton step was not good; reduce the radius */ 698 tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 699 } 700 break; 701 702 case BNK_UPDATE_REDUCTION: 703 if (stepType == BNK_NEWTON) { 704 if (prered < 0.0) { 705 /* The predicted reduction has the wrong sign. This cannot 706 happen in infinite precision arithmetic. Step should 707 be rejected! */ 708 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 709 } 710 else { 711 if (PetscIsInfOrNanReal(actred)) { 712 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 713 } else { 714 if ((PetscAbsScalar(actred) <= bnk->epsilon) && 715 (PetscAbsScalar(prered) <= bnk->epsilon)) { 716 kappa = 1.0; 717 } 718 else { 719 kappa = actred / prered; 720 } 721 722 /* Accept or reject the step and update radius */ 723 if (kappa < bnk->eta1) { 724 /* Reject the step */ 725 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 726 } 727 else { 728 /* Accept the step */ 729 if (kappa < bnk->eta2) { 730 /* Marginal bad step */ 731 tao->trust = bnk->alpha2 * PetscMin(tao->trust, bnk->dnorm); 732 } 733 else if (kappa < bnk->eta3) { 734 /* Reasonable step */ 735 tao->trust = bnk->alpha3 * tao->trust; 736 } 737 else if (kappa < bnk->eta4) { 738 /* Good step */ 739 tao->trust = PetscMax(bnk->alpha4 * bnk->dnorm, tao->trust); 740 } 741 else { 742 /* Very good step */ 743 tao->trust = PetscMax(bnk->alpha5 * bnk->dnorm, tao->trust); 744 } 745 *accept = PETSC_TRUE; 746 } 747 } 748 } 749 } else { 750 /* Newton step was not good; reduce the radius */ 751 tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 752 } 753 break; 754 755 default: 756 if (stepType == BNK_NEWTON) { 757 if (prered < 0.0) { 758 /* The predicted reduction has the wrong sign. This cannot */ 759 /* happen in infinite precision arithmetic. Step should */ 760 /* be rejected! */ 761 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 762 } else { 763 if (PetscIsInfOrNanReal(actred)) { 764 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 765 } else { 766 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 767 kappa = 1.0; 768 } else { 769 kappa = actred / prered; 770 } 771 772 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 773 tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 774 tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 775 tau_min = PetscMin(tau_1, tau_2); 776 tau_max = PetscMax(tau_1, tau_2); 777 778 if (kappa >= 1.0 - bnk->mu1) { 779 /* Great agreement */ 780 *accept = PETSC_TRUE; 781 if (tau_max < 1.0) { 782 tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 783 } else if (tau_max > bnk->gamma4) { 784 tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 785 } else { 786 tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 787 } 788 } else if (kappa >= 1.0 - bnk->mu2) { 789 /* Good agreement */ 790 *accept = PETSC_TRUE; 791 if (tau_max < bnk->gamma2) { 792 tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 793 } else if (tau_max > bnk->gamma3) { 794 tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 795 } else if (tau_max < 1.0) { 796 tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 797 } else { 798 tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 799 } 800 } else { 801 /* Not good agreement */ 802 if (tau_min > 1.0) { 803 tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 804 } else if (tau_max < bnk->gamma1) { 805 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 806 } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 807 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 808 } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 809 tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 810 } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 811 tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 812 } else { 813 tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 814 } 815 } 816 } 817 } 818 } else { 819 /* Newton step was not good; reduce the radius */ 820 tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 821 } 822 /* The radius may have been increased; modify if it is too large */ 823 tao->trust = PetscMin(tao->trust, bnk->max_radius); 824 } 825 tao->trust = PetscMax(tao->trust, bnk->min_radius); 826 PetscFunctionReturn(0); 827 } 828 829 /* ---------------------------------------------------------- */ 830 831 static PetscErrorCode TaoSetUp_BNK(Tao tao) 832 { 833 TAO_BNK *bnk = (TAO_BNK *)tao->data; 834 PetscErrorCode ierr; 835 836 PetscFunctionBegin; 837 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 838 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 839 if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 840 if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 841 if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 842 if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);} 843 if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);} 844 if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);} 845 if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);} 846 bnk->Diag = 0; 847 bnk->M = 0; 848 bnk->H_inactive = 0; 849 bnk->Hpre_inactive = 0; 850 PetscFunctionReturn(0); 851 } 852 853 /*------------------------------------------------------------*/ 854 855 static PetscErrorCode TaoDestroy_BNK(Tao tao) 856 { 857 TAO_BNK *bnk = (TAO_BNK *)tao->data; 858 PetscErrorCode ierr; 859 860 PetscFunctionBegin; 861 if (tao->setupcalled) { 862 ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 863 ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 864 ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 865 ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr); 866 ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr); 867 ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr); 868 ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr); 869 } 870 ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 871 ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 872 ierr = PetscFree(tao->data);CHKERRQ(ierr); 873 PetscFunctionReturn(0); 874 } 875 876 /*------------------------------------------------------------*/ 877 878 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 879 { 880 TAO_BNK *bnk = (TAO_BNK *)tao->data; 881 PetscErrorCode ierr; 882 883 PetscFunctionBegin; 884 ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 885 ierr = PetscOptionsEList("-tao_BNK_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr); 886 ierr = PetscOptionsEList("-tao_BNK_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[bnk->bfgs_scale_type], &bnk->bfgs_scale_type, 0);CHKERRQ(ierr); 887 ierr = PetscOptionsEList("-tao_BNK_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, 0);CHKERRQ(ierr); 888 ierr = PetscOptionsEList("-tao_BNK_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);CHKERRQ(ierr); 889 ierr = PetscOptionsReal("-tao_BNK_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 890 ierr = PetscOptionsReal("-tao_BNK_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 891 ierr = PetscOptionsReal("-tao_BNK_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 892 ierr = PetscOptionsReal("-tao_BNK_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 893 ierr = PetscOptionsReal("-tao_BNK_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 894 ierr = PetscOptionsReal("-tao_BNK_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 895 ierr = PetscOptionsReal("-tao_BNK_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 896 ierr = PetscOptionsReal("-tao_BNK_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 897 ierr = PetscOptionsReal("-tao_BNK_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 898 ierr = PetscOptionsReal("-tao_BNK_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 899 ierr = PetscOptionsReal("-tao_BNK_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 900 ierr = PetscOptionsReal("-tao_BNK_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 901 ierr = PetscOptionsReal("-tao_BNK_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 902 ierr = PetscOptionsReal("-tao_BNK_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 903 ierr = PetscOptionsReal("-tao_BNK_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 904 ierr = PetscOptionsReal("-tao_BNK_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 905 ierr = PetscOptionsReal("-tao_BNK_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 906 ierr = PetscOptionsReal("-tao_BNK_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 907 ierr = PetscOptionsReal("-tao_BNK_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 908 ierr = PetscOptionsReal("-tao_BNK_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 909 ierr = PetscOptionsReal("-tao_BNK_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 910 ierr = PetscOptionsReal("-tao_BNK_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 911 ierr = PetscOptionsReal("-tao_BNK_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 912 ierr = PetscOptionsReal("-tao_BNK_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 913 ierr = PetscOptionsReal("-tao_BNK_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 914 ierr = PetscOptionsReal("-tao_BNK_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 915 ierr = PetscOptionsReal("-tao_BNK_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 916 ierr = PetscOptionsReal("-tao_BNK_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 917 ierr = PetscOptionsReal("-tao_BNK_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 918 ierr = PetscOptionsReal("-tao_BNK_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 919 ierr = PetscOptionsReal("-tao_BNK_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 920 ierr = PetscOptionsReal("-tao_BNK_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 921 ierr = PetscOptionsReal("-tao_BNK_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 922 ierr = PetscOptionsReal("-tao_BNK_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 923 ierr = PetscOptionsReal("-tao_BNK_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 924 ierr = PetscOptionsReal("-tao_BNK_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 925 ierr = PetscOptionsReal("-tao_BNK_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 926 ierr = PetscOptionsReal("-tao_BNK_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 927 ierr = PetscOptionsReal("-tao_BNK_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 928 ierr = PetscOptionsReal("-tao_BNK_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 929 ierr = PetscOptionsReal("-tao_BNK_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 930 ierr = PetscOptionsReal("-tao_BNK_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 931 ierr = PetscOptionsReal("-tao_BNK_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 932 ierr = PetscOptionsReal("-tao_BNK_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 933 ierr = PetscOptionsReal("-tao_BNK_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 934 ierr = PetscOptionsTail();CHKERRQ(ierr); 935 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 936 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 937 PetscFunctionReturn(0); 938 } 939 940 /*------------------------------------------------------------*/ 941 942 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 943 { 944 TAO_BNK *bnk = (TAO_BNK *)tao->data; 945 PetscInt nrejects; 946 PetscBool isascii; 947 PetscErrorCode ierr; 948 949 PetscFunctionBegin; 950 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 951 if (isascii) { 952 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 953 if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 954 ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 955 ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 956 } 957 ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 958 ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 959 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 960 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 961 ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 962 ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 963 ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 964 ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 965 ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 966 ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 967 ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 968 ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 969 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 970 } 971 PetscFunctionReturn(0); 972 } 973 974 /* ---------------------------------------------------------- */ 975 976 /*MC 977 TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 978 At each iteration, the BNK methods solve the symmetric 979 system of equations to obtain the step diretion dk: 980 Hk dk = -gk 981 for free variables only. The step can be globalized either through 982 trust-region methods, or a line search, or a heuristic mixture of both. 983 984 Options Database Keys: 985 + -tao_BNK_pc_type - "none","ahess","bfgs","petsc" 986 . -tao_BNK_bfgs_scale_type - "ahess","phess","bfgs" 987 . -tao_BNK_init_type - "constant","direction","interpolation" 988 . -tao_BNK_update_type - "step","direction","interpolation" 989 . -tao_BNK_sval - perturbation starting value 990 . -tao_BNK_imin - minimum initial perturbation 991 . -tao_BNK_imax - maximum initial perturbation 992 . -tao_BNK_pmin - minimum perturbation 993 . -tao_BNK_pmax - maximum perturbation 994 . -tao_BNK_pgfac - growth factor 995 . -tao_BNK_psfac - shrink factor 996 . -tao_BNK_imfac - initial merit factor 997 . -tao_BNK_pmgfac - merit growth factor 998 . -tao_BNK_pmsfac - merit shrink factor 999 . -tao_BNK_eta1 - poor steplength; reduce radius 1000 . -tao_BNK_eta2 - reasonable steplength; leave radius 1001 . -tao_BNK_eta3 - good steplength; increase readius 1002 . -tao_BNK_eta4 - excellent steplength; greatly increase radius 1003 . -tao_BNK_alpha1 - alpha1 reduction 1004 . -tao_BNK_alpha2 - alpha2 reduction 1005 . -tao_BNK_alpha3 - alpha3 reduction 1006 . -tao_BNK_alpha4 - alpha4 reduction 1007 . -tao_BNK_alpha - alpha5 reduction 1008 . -tao_BNK_mu1 - mu1 interpolation update 1009 . -tao_BNK_mu2 - mu2 interpolation update 1010 . -tao_BNK_gamma1 - gamma1 interpolation update 1011 . -tao_BNK_gamma2 - gamma2 interpolation update 1012 . -tao_BNK_gamma3 - gamma3 interpolation update 1013 . -tao_BNK_gamma4 - gamma4 interpolation update 1014 . -tao_BNK_theta - theta interpolation update 1015 . -tao_BNK_omega1 - omega1 step update 1016 . -tao_BNK_omega2 - omega2 step update 1017 . -tao_BNK_omega3 - omega3 step update 1018 . -tao_BNK_omega4 - omega4 step update 1019 . -tao_BNK_omega5 - omega5 step update 1020 . -tao_BNK_mu1_i - mu1 interpolation init factor 1021 . -tao_BNK_mu2_i - mu2 interpolation init factor 1022 . -tao_BNK_gamma1_i - gamma1 interpolation init factor 1023 . -tao_BNK_gamma2_i - gamma2 interpolation init factor 1024 . -tao_BNK_gamma3_i - gamma3 interpolation init factor 1025 . -tao_BNK_gamma4_i - gamma4 interpolation init factor 1026 - -tao_BNK_theta_i - theta interpolation init factor 1027 1028 Level: beginner 1029 M*/ 1030 1031 PetscErrorCode TaoCreate_BNK(Tao tao) 1032 { 1033 TAO_BNK *bnk; 1034 const char *morethuente_type = TAOLINESEARCHMT; 1035 PetscErrorCode ierr; 1036 1037 PetscFunctionBegin; 1038 ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 1039 1040 tao->ops->setup = TaoSetUp_BNK; 1041 tao->ops->view = TaoView_BNK; 1042 tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1043 tao->ops->destroy = TaoDestroy_BNK; 1044 1045 /* Override default settings (unless already changed) */ 1046 if (!tao->max_it_changed) tao->max_it = 50; 1047 if (!tao->trust0_changed) tao->trust0 = 100.0; 1048 1049 tao->data = (void*)bnk; 1050 1051 /* Hessian shifting parameters */ 1052 bnk->sval = 0.0; 1053 bnk->imin = 1.0e-4; 1054 bnk->imax = 1.0e+2; 1055 bnk->imfac = 1.0e-1; 1056 1057 bnk->pmin = 1.0e-12; 1058 bnk->pmax = 1.0e+2; 1059 bnk->pgfac = 1.0e+1; 1060 bnk->psfac = 4.0e-1; 1061 bnk->pmgfac = 1.0e-1; 1062 bnk->pmsfac = 1.0e-1; 1063 1064 /* Default values for trust-region radius update based on steplength */ 1065 bnk->nu1 = 0.25; 1066 bnk->nu2 = 0.50; 1067 bnk->nu3 = 1.00; 1068 bnk->nu4 = 1.25; 1069 1070 bnk->omega1 = 0.25; 1071 bnk->omega2 = 0.50; 1072 bnk->omega3 = 1.00; 1073 bnk->omega4 = 2.00; 1074 bnk->omega5 = 4.00; 1075 1076 /* Default values for trust-region radius update based on reduction */ 1077 bnk->eta1 = 1.0e-4; 1078 bnk->eta2 = 0.25; 1079 bnk->eta3 = 0.50; 1080 bnk->eta4 = 0.90; 1081 1082 bnk->alpha1 = 0.25; 1083 bnk->alpha2 = 0.50; 1084 bnk->alpha3 = 1.00; 1085 bnk->alpha4 = 2.00; 1086 bnk->alpha5 = 4.00; 1087 1088 /* Default values for trust-region radius update based on interpolation */ 1089 bnk->mu1 = 0.10; 1090 bnk->mu2 = 0.50; 1091 1092 bnk->gamma1 = 0.25; 1093 bnk->gamma2 = 0.50; 1094 bnk->gamma3 = 2.00; 1095 bnk->gamma4 = 4.00; 1096 1097 bnk->theta = 0.05; 1098 1099 /* Default values for trust region initialization based on interpolation */ 1100 bnk->mu1_i = 0.35; 1101 bnk->mu2_i = 0.50; 1102 1103 bnk->gamma1_i = 0.0625; 1104 bnk->gamma2_i = 0.5; 1105 bnk->gamma3_i = 2.0; 1106 bnk->gamma4_i = 5.0; 1107 1108 bnk->theta_i = 0.25; 1109 1110 /* Remaining parameters */ 1111 bnk->min_radius = 1.0e-10; 1112 bnk->max_radius = 1.0e10; 1113 bnk->epsilon = 1.0e-6; 1114 1115 bnk->pc_type = BNK_PC_BFGS; 1116 bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 1117 bnk->init_type = BNK_INIT_INTERPOLATION; 1118 bnk->update_type = BNK_UPDATE_INTERPOLATION; 1119 1120 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1121 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 1122 ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1123 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1124 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1125 1126 /* Set linear solver to default for symmetric matrices */ 1127 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1128 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1129 ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 1130 ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 1131 PetscFunctionReturn(0); 1132 } 1133