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