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