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