1 #include <../src/tao/unconstrained/impls/ntl/ntlimpl.h> 2 3 #include <petscksp.h> 4 5 #define NTL_INIT_CONSTANT 0 6 #define NTL_INIT_DIRECTION 1 7 #define NTL_INIT_INTERPOLATION 2 8 #define NTL_INIT_TYPES 3 9 10 #define NTL_UPDATE_REDUCTION 0 11 #define NTL_UPDATE_INTERPOLATION 1 12 #define NTL_UPDATE_TYPES 2 13 14 static const char *NTL_INIT[64] = {"constant","direction","interpolation"}; 15 16 static const char *NTL_UPDATE[64] = {"reduction","interpolation"}; 17 18 /* Implements Newton's Method with a trust-region, line-search approach for 19 solving unconstrained minimization problems. A More'-Thuente line search 20 is used to guarantee that the bfgs preconditioner remains positive 21 definite. */ 22 23 #define NTL_NEWTON 0 24 #define NTL_BFGS 1 25 #define NTL_SCALED_GRADIENT 2 26 #define NTL_GRADIENT 3 27 28 static PetscErrorCode TaoSolve_NTL(Tao tao) 29 { 30 TAO_NTL *tl = (TAO_NTL *)tao->data; 31 KSPType ksp_type; 32 PetscBool is_nash,is_stcg,is_gltr,is_bfgs,is_jacobi; 33 KSPConvergedReason ksp_reason; 34 PC pc; 35 TaoLineSearchConvergedReason ls_reason; 36 37 PetscReal fmin, ftrial, prered, actred, kappa, sigma; 38 PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 39 PetscReal f, fold, gdx, gnorm; 40 PetscReal step = 1.0; 41 42 PetscReal norm_d = 0.0; 43 PetscErrorCode ierr; 44 PetscInt stepType; 45 PetscInt its; 46 47 PetscInt bfgsUpdates = 0; 48 PetscInt needH; 49 50 PetscInt i_max = 5; 51 PetscInt j_max = 1; 52 PetscInt i, j, n, N; 53 54 PetscInt tr_reject; 55 56 PetscFunctionBegin; 57 if (tao->XL || tao->XU || tao->ops->computebounds) { 58 ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by ntl algorithm\n");CHKERRQ(ierr); 59 } 60 61 ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 62 ierr = PetscStrcmp(ksp_type,KSPCGNASH,&is_nash);CHKERRQ(ierr); 63 ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&is_stcg);CHKERRQ(ierr); 64 ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&is_gltr);CHKERRQ(ierr); 65 if (!is_nash && !is_stcg && !is_gltr) { 66 SETERRQ(PETSC_COMM_SELF,1,"TAO_NTR requires nash, stcg, or gltr for the KSP"); 67 } 68 69 /* Initialize the radius and modify if it is too large or small */ 70 tao->trust = tao->trust0; 71 tao->trust = PetscMax(tao->trust, tl->min_radius); 72 tao->trust = PetscMin(tao->trust, tl->max_radius); 73 74 /* Allocate the vectors needed for the BFGS approximation */ 75 ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 76 ierr = PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);CHKERRQ(ierr); 77 ierr = PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);CHKERRQ(ierr); 78 if (is_bfgs) { 79 tl->bfgs_pre = pc; 80 ierr = PCLMVMGetMatLMVM(tl->bfgs_pre, &tl->M);CHKERRQ(ierr); 81 ierr = VecGetLocalSize(tao->solution, &n);CHKERRQ(ierr); 82 ierr = VecGetSize(tao->solution, &N);CHKERRQ(ierr); 83 ierr = MatSetSizes(tl->M, n, n, N, N);CHKERRQ(ierr); 84 ierr = MatLMVMAllocate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 85 } else if (is_jacobi) { 86 ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 87 } 88 89 /* Check convergence criteria */ 90 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 91 ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 92 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 93 needH = 1; 94 95 tao->reason = TAO_CONTINUE_ITERATING; 96 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 97 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 98 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 99 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 100 101 /* Initialize trust-region radius */ 102 switch(tl->init_type) { 103 case NTL_INIT_CONSTANT: 104 /* Use the initial radius specified */ 105 break; 106 107 case NTL_INIT_INTERPOLATION: 108 /* Use the initial radius specified */ 109 max_radius = 0.0; 110 111 for (j = 0; j < j_max; ++j) { 112 fmin = f; 113 sigma = 0.0; 114 115 if (needH) { 116 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 117 needH = 0; 118 } 119 120 for (i = 0; i < i_max; ++i) { 121 ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 122 ierr = VecAXPY(tl->W, -tao->trust/gnorm, tao->gradient);CHKERRQ(ierr); 123 124 ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 125 if (PetscIsInfOrNanReal(ftrial)) { 126 tau = tl->gamma1_i; 127 } else { 128 if (ftrial < fmin) { 129 fmin = ftrial; 130 sigma = -tao->trust / gnorm; 131 } 132 133 ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 134 ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr); 135 136 prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 137 actred = f - ftrial; 138 if ((PetscAbsScalar(actred) <= tl->epsilon) && (PetscAbsScalar(prered) <= tl->epsilon)) { 139 kappa = 1.0; 140 } else { 141 kappa = actred / prered; 142 } 143 144 tau_1 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust + (1.0 - tl->theta_i) * prered - actred); 145 tau_2 = tl->theta_i * gnorm * tao->trust / (tl->theta_i * gnorm * tao->trust - (1.0 + tl->theta_i) * prered + actred); 146 tau_min = PetscMin(tau_1, tau_2); 147 tau_max = PetscMax(tau_1, tau_2); 148 149 if (PetscAbsScalar(kappa - 1.0) <= tl->mu1_i) { 150 /* Great agreement */ 151 max_radius = PetscMax(max_radius, tao->trust); 152 153 if (tau_max < 1.0) { 154 tau = tl->gamma3_i; 155 } else if (tau_max > tl->gamma4_i) { 156 tau = tl->gamma4_i; 157 } else if (tau_1 >= 1.0 && tau_1 <= tl->gamma4_i && tau_2 < 1.0) { 158 tau = tau_1; 159 } else if (tau_2 >= 1.0 && tau_2 <= tl->gamma4_i && tau_1 < 1.0) { 160 tau = tau_2; 161 } else { 162 tau = tau_max; 163 } 164 } else if (PetscAbsScalar(kappa - 1.0) <= tl->mu2_i) { 165 /* Good agreement */ 166 max_radius = PetscMax(max_radius, tao->trust); 167 168 if (tau_max < tl->gamma2_i) { 169 tau = tl->gamma2_i; 170 } else if (tau_max > tl->gamma3_i) { 171 tau = tl->gamma3_i; 172 } else { 173 tau = tau_max; 174 } 175 } else { 176 /* Not good agreement */ 177 if (tau_min > 1.0) { 178 tau = tl->gamma2_i; 179 } else if (tau_max < tl->gamma1_i) { 180 tau = tl->gamma1_i; 181 } else if ((tau_min < tl->gamma1_i) && (tau_max >= 1.0)) { 182 tau = tl->gamma1_i; 183 } else if ((tau_1 >= tl->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1_i) || (tau_2 >= 1.0))) { 184 tau = tau_1; 185 } else if ((tau_2 >= tl->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1_i) || (tau_2 >= 1.0))) { 186 tau = tau_2; 187 } else { 188 tau = tau_max; 189 } 190 } 191 } 192 tao->trust = tau * tao->trust; 193 } 194 195 if (fmin < f) { 196 f = fmin; 197 ierr = VecAXPY(tao->solution, sigma, tao->gradient);CHKERRQ(ierr); 198 ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 199 200 ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 201 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 202 needH = 1; 203 204 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 205 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 206 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 207 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 208 } 209 } 210 tao->trust = PetscMax(tao->trust, max_radius); 211 212 /* Modify the radius if it is too large or small */ 213 tao->trust = PetscMax(tao->trust, tl->min_radius); 214 tao->trust = PetscMin(tao->trust, tl->max_radius); 215 break; 216 217 default: 218 /* Norm of the first direction will initialize radius */ 219 tao->trust = 0.0; 220 break; 221 } 222 223 /* Set counter for gradient/reset steps */ 224 tl->ntrust = 0; 225 tl->newt = 0; 226 tl->bfgs = 0; 227 tl->grad = 0; 228 229 /* Have not converged; continue with Newton method */ 230 while (tao->reason == TAO_CONTINUE_ITERATING) { 231 ++tao->niter; 232 tao->ksp_its=0; 233 /* Compute the Hessian */ 234 if (needH) { 235 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 236 } 237 238 if (tl->bfgs_pre) { 239 /* Update the limited memory preconditioner */ 240 ierr = MatLMVMUpdate(tl->M,tao->solution, tao->gradient);CHKERRQ(ierr); 241 ++bfgsUpdates; 242 } 243 ierr = KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre);CHKERRQ(ierr); 244 /* Solve the Newton system of equations */ 245 ierr = KSPCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 246 ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 247 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 248 tao->ksp_its+=its; 249 tao->ksp_tot_its+=its; 250 ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 251 252 if (0.0 == tao->trust) { 253 /* Radius was uninitialized; use the norm of the direction */ 254 if (norm_d > 0.0) { 255 tao->trust = norm_d; 256 257 /* Modify the radius if it is too large or small */ 258 tao->trust = PetscMax(tao->trust, tl->min_radius); 259 tao->trust = PetscMin(tao->trust, tl->max_radius); 260 } else { 261 /* The direction was bad; set radius to default value and re-solve 262 the trust-region subproblem to get a direction */ 263 tao->trust = tao->trust0; 264 265 /* Modify the radius if it is too large or small */ 266 tao->trust = PetscMax(tao->trust, tl->min_radius); 267 tao->trust = PetscMin(tao->trust, tl->max_radius); 268 269 ierr = KSPCGSetRadius(tao->ksp,tl->max_radius);CHKERRQ(ierr); 270 ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 271 ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 272 tao->ksp_its+=its; 273 tao->ksp_tot_its+=its; 274 ierr = KSPCGGetNormD(tao->ksp, &norm_d);CHKERRQ(ierr); 275 276 if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 277 } 278 } 279 280 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 281 ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 282 if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (tl->bfgs_pre)) { 283 /* Preconditioner is numerically indefinite; reset the 284 approximate if using BFGS preconditioning. */ 285 ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 286 ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 287 bfgsUpdates = 1; 288 } 289 290 /* Check trust-region reduction conditions */ 291 tr_reject = 0; 292 if (NTL_UPDATE_REDUCTION == tl->update_type) { 293 /* Get predicted reduction */ 294 ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 295 if (prered >= 0.0) { 296 /* The predicted reduction has the wrong sign. This cannot 297 happen in infinite precision arithmetic. Step should 298 be rejected! */ 299 tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 300 tr_reject = 1; 301 } else { 302 /* Compute trial step and function value */ 303 ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 304 ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 305 ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 306 307 if (PetscIsInfOrNanReal(ftrial)) { 308 tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 309 tr_reject = 1; 310 } else { 311 /* Compute and actual reduction */ 312 actred = f - ftrial; 313 prered = -prered; 314 if ((PetscAbsScalar(actred) <= tl->epsilon) && 315 (PetscAbsScalar(prered) <= tl->epsilon)) { 316 kappa = 1.0; 317 } else { 318 kappa = actred / prered; 319 } 320 321 /* Accept of reject the step and update radius */ 322 if (kappa < tl->eta1) { 323 /* Reject the step */ 324 tao->trust = tl->alpha1 * PetscMin(tao->trust, norm_d); 325 tr_reject = 1; 326 } else { 327 /* Accept the step */ 328 if (kappa < tl->eta2) { 329 /* Marginal bad step */ 330 tao->trust = tl->alpha2 * PetscMin(tao->trust, norm_d); 331 } else if (kappa < tl->eta3) { 332 /* Reasonable step */ 333 tao->trust = tl->alpha3 * tao->trust; 334 } else if (kappa < tl->eta4) { 335 /* Good step */ 336 tao->trust = PetscMax(tl->alpha4 * norm_d, tao->trust); 337 } else { 338 /* Very good step */ 339 tao->trust = PetscMax(tl->alpha5 * norm_d, tao->trust); 340 } 341 } 342 } 343 } 344 } else { 345 /* Get predicted reduction */ 346 ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 347 if (prered >= 0.0) { 348 /* The predicted reduction has the wrong sign. This cannot 349 happen in infinite precision arithmetic. Step should 350 be rejected! */ 351 tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 352 tr_reject = 1; 353 } else { 354 ierr = VecCopy(tao->solution, tl->W);CHKERRQ(ierr); 355 ierr = VecAXPY(tl->W, 1.0, tao->stepdirection);CHKERRQ(ierr); 356 ierr = TaoComputeObjective(tao, tl->W, &ftrial);CHKERRQ(ierr); 357 if (PetscIsInfOrNanReal(ftrial)) { 358 tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 359 tr_reject = 1; 360 } else { 361 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 362 363 actred = f - ftrial; 364 prered = -prered; 365 if ((PetscAbsScalar(actred) <= tl->epsilon) && 366 (PetscAbsScalar(prered) <= tl->epsilon)) { 367 kappa = 1.0; 368 } else { 369 kappa = actred / prered; 370 } 371 372 tau_1 = tl->theta * gdx / (tl->theta * gdx - (1.0 - tl->theta) * prered + actred); 373 tau_2 = tl->theta * gdx / (tl->theta * gdx + (1.0 + tl->theta) * prered - actred); 374 tau_min = PetscMin(tau_1, tau_2); 375 tau_max = PetscMax(tau_1, tau_2); 376 377 if (kappa >= 1.0 - tl->mu1) { 378 /* Great agreement; accept step and update radius */ 379 if (tau_max < 1.0) { 380 tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 381 } else if (tau_max > tl->gamma4) { 382 tao->trust = PetscMax(tao->trust, tl->gamma4 * norm_d); 383 } else { 384 tao->trust = PetscMax(tao->trust, tau_max * norm_d); 385 } 386 } else if (kappa >= 1.0 - tl->mu2) { 387 /* Good agreement */ 388 389 if (tau_max < tl->gamma2) { 390 tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 391 } else if (tau_max > tl->gamma3) { 392 tao->trust = PetscMax(tao->trust, tl->gamma3 * norm_d); 393 } else if (tau_max < 1.0) { 394 tao->trust = tau_max * PetscMin(tao->trust, norm_d); 395 } else { 396 tao->trust = PetscMax(tao->trust, tau_max * norm_d); 397 } 398 } else { 399 /* Not good agreement */ 400 if (tau_min > 1.0) { 401 tao->trust = tl->gamma2 * PetscMin(tao->trust, norm_d); 402 } else if (tau_max < tl->gamma1) { 403 tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 404 } else if ((tau_min < tl->gamma1) && (tau_max >= 1.0)) { 405 tao->trust = tl->gamma1 * PetscMin(tao->trust, norm_d); 406 } else if ((tau_1 >= tl->gamma1) && (tau_1 < 1.0) && ((tau_2 < tl->gamma1) || (tau_2 >= 1.0))) { 407 tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 408 } else if ((tau_2 >= tl->gamma1) && (tau_2 < 1.0) && ((tau_1 < tl->gamma1) || (tau_2 >= 1.0))) { 409 tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 410 } else { 411 tao->trust = tau_max * PetscMin(tao->trust, norm_d); 412 } 413 tr_reject = 1; 414 } 415 } 416 } 417 } 418 419 if (tr_reject) { 420 /* The trust-region constraints rejected the step. Apply a linesearch. 421 Check for descent direction. */ 422 ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 423 if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 424 /* Newton step is not descent or direction produced Inf or NaN */ 425 426 if (!tl->bfgs_pre) { 427 /* We don't have the bfgs matrix around and updated 428 Must use gradient direction in this case */ 429 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 430 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 431 ++tl->grad; 432 stepType = NTL_GRADIENT; 433 } else { 434 /* Attempt to use the BFGS direction */ 435 ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 436 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 437 438 /* Check for success (descent direction) */ 439 ierr = VecDot(tao->stepdirection, tao->gradient, &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 ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 448 ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 449 ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 450 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 451 452 bfgsUpdates = 1; 453 ++tl->grad; 454 stepType = NTL_GRADIENT; 455 } else { 456 ierr = MatLMVMGetUpdateCount(tl->M, &bfgsUpdates);CHKERRQ(ierr); 457 if (1 == bfgsUpdates) { 458 /* The first BFGS direction is always the scaled gradient */ 459 ++tl->grad; 460 stepType = NTL_GRADIENT; 461 } else { 462 ++tl->bfgs; 463 stepType = NTL_BFGS; 464 } 465 } 466 } 467 } else { 468 /* Computed Newton step is descent */ 469 ++tl->newt; 470 stepType = NTL_NEWTON; 471 } 472 473 /* Perform the linesearch */ 474 fold = f; 475 ierr = VecCopy(tao->solution, tl->Xold);CHKERRQ(ierr); 476 ierr = VecCopy(tao->gradient, tl->Gold);CHKERRQ(ierr); 477 478 step = 1.0; 479 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 480 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 481 482 while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NTL_GRADIENT) { /* Linesearch failed */ 483 /* Linesearch failed */ 484 f = fold; 485 ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 486 ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 487 488 switch(stepType) { 489 case NTL_NEWTON: 490 /* Failed to obtain acceptable iterate with Newton step */ 491 492 if (tl->bfgs_pre) { 493 /* We don't have the bfgs matrix around and being updated 494 Must use gradient direction in this case */ 495 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 496 ++tl->grad; 497 stepType = NTL_GRADIENT; 498 } else { 499 /* Attempt to use the BFGS direction */ 500 ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 501 502 503 /* Check for success (descent direction) */ 504 ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 505 if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 506 /* BFGS direction is not descent or direction produced 507 not a number. We can assert bfgsUpdates > 1 in this case 508 Use steepest descent direction (scaled) */ 509 ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 510 ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 511 ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 512 513 bfgsUpdates = 1; 514 ++tl->grad; 515 stepType = NTL_GRADIENT; 516 } else { 517 ierr = MatLMVMGetUpdateCount(tl->M, &bfgsUpdates);CHKERRQ(ierr); 518 if (1 == bfgsUpdates) { 519 /* The first BFGS direction is always the scaled gradient */ 520 ++tl->grad; 521 stepType = NTL_GRADIENT; 522 } else { 523 ++tl->bfgs; 524 stepType = NTL_BFGS; 525 } 526 } 527 } 528 break; 529 530 case NTL_BFGS: 531 /* Can only enter if pc_type == NTL_PC_BFGS 532 Failed to obtain acceptable iterate with BFGS step 533 Attempt to use the scaled gradient direction */ 534 ierr = MatLMVMReset(tl->M, PETSC_FALSE);CHKERRQ(ierr); 535 ierr = MatLMVMUpdate(tl->M, tao->solution, tao->gradient);CHKERRQ(ierr); 536 ierr = MatSolve(tl->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 537 538 bfgsUpdates = 1; 539 ++tl->grad; 540 stepType = NTL_GRADIENT; 541 break; 542 } 543 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 544 545 /* This may be incorrect; linesearch has values for stepmax and stepmin 546 that should be reset. */ 547 step = 1.0; 548 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 549 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 550 } 551 552 if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 553 /* Failed to find an improving point */ 554 f = fold; 555 ierr = VecCopy(tl->Xold, tao->solution);CHKERRQ(ierr); 556 ierr = VecCopy(tl->Gold, tao->gradient);CHKERRQ(ierr); 557 tao->trust = 0.0; 558 step = 0.0; 559 tao->reason = TAO_DIVERGED_LS_FAILURE; 560 break; 561 } else if (stepType == NTL_NEWTON) { 562 if (step < tl->nu1) { 563 /* Very bad step taken; reduce radius */ 564 tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 565 } else if (step < tl->nu2) { 566 /* Reasonably bad step taken; reduce radius */ 567 tao->trust = tl->omega2 * PetscMin(norm_d, tao->trust); 568 } else if (step < tl->nu3) { 569 /* Reasonable step was taken; leave radius alone */ 570 if (tl->omega3 < 1.0) { 571 tao->trust = tl->omega3 * PetscMin(norm_d, tao->trust); 572 } else if (tl->omega3 > 1.0) { 573 tao->trust = PetscMax(tl->omega3 * norm_d, tao->trust); 574 } 575 } else if (step < tl->nu4) { 576 /* Full step taken; increase the radius */ 577 tao->trust = PetscMax(tl->omega4 * norm_d, tao->trust); 578 } else { 579 /* More than full step taken; increase the radius */ 580 tao->trust = PetscMax(tl->omega5 * norm_d, tao->trust); 581 } 582 } else { 583 /* Newton step was not good; reduce the radius */ 584 tao->trust = tl->omega1 * PetscMin(norm_d, tao->trust); 585 } 586 } else { 587 /* Trust-region step is accepted */ 588 ierr = VecCopy(tl->W, tao->solution);CHKERRQ(ierr); 589 f = ftrial; 590 ierr = TaoComputeGradient(tao, tao->solution, tao->gradient);CHKERRQ(ierr); 591 ++tl->ntrust; 592 } 593 594 /* The radius may have been increased; modify if it is too large */ 595 tao->trust = PetscMin(tao->trust, tl->max_radius); 596 597 /* Check for converged */ 598 ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 599 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 600 needH = 1; 601 602 ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 603 ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step);CHKERRQ(ierr); 604 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 605 } 606 PetscFunctionReturn(0); 607 } 608 609 /* ---------------------------------------------------------- */ 610 static PetscErrorCode TaoSetUp_NTL(Tao tao) 611 { 612 TAO_NTL *tl = (TAO_NTL *)tao->data; 613 PetscErrorCode ierr; 614 615 PetscFunctionBegin; 616 if (!tao->gradient) {ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); } 617 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);} 618 if (!tl->W) { ierr = VecDuplicate(tao->solution, &tl->W);CHKERRQ(ierr);} 619 if (!tl->Xold) { ierr = VecDuplicate(tao->solution, &tl->Xold);CHKERRQ(ierr);} 620 if (!tl->Gold) { ierr = VecDuplicate(tao->solution, &tl->Gold);CHKERRQ(ierr);} 621 tl->bfgs_pre = 0; 622 tl->M = 0; 623 PetscFunctionReturn(0); 624 } 625 626 /*------------------------------------------------------------*/ 627 static PetscErrorCode TaoDestroy_NTL(Tao tao) 628 { 629 TAO_NTL *tl = (TAO_NTL *)tao->data; 630 PetscErrorCode ierr; 631 632 PetscFunctionBegin; 633 if (tao->setupcalled) { 634 ierr = VecDestroy(&tl->W);CHKERRQ(ierr); 635 ierr = VecDestroy(&tl->Xold);CHKERRQ(ierr); 636 ierr = VecDestroy(&tl->Gold);CHKERRQ(ierr); 637 } 638 ierr = PetscFree(tao->data);CHKERRQ(ierr); 639 PetscFunctionReturn(0); 640 } 641 642 /*------------------------------------------------------------*/ 643 static PetscErrorCode TaoSetFromOptions_NTL(PetscOptionItems *PetscOptionsObject,Tao tao) 644 { 645 TAO_NTL *tl = (TAO_NTL *)tao->data; 646 PetscErrorCode ierr; 647 648 PetscFunctionBegin; 649 ierr = PetscOptionsHead(PetscOptionsObject,"Newton trust region with line search method for unconstrained optimization");CHKERRQ(ierr); 650 ierr = PetscOptionsEList("-tao_ntl_init_type", "radius initialization type", "", NTL_INIT, NTL_INIT_TYPES, NTL_INIT[tl->init_type], &tl->init_type,NULL);CHKERRQ(ierr); 651 ierr = PetscOptionsEList("-tao_ntl_update_type", "radius update type", "", NTL_UPDATE, NTL_UPDATE_TYPES, NTL_UPDATE[tl->update_type], &tl->update_type,NULL);CHKERRQ(ierr); 652 ierr = PetscOptionsReal("-tao_ntl_eta1", "poor steplength; reduce radius", "", tl->eta1, &tl->eta1,NULL);CHKERRQ(ierr); 653 ierr = PetscOptionsReal("-tao_ntl_eta2", "reasonable steplength; leave radius alone", "", tl->eta2, &tl->eta2,NULL);CHKERRQ(ierr); 654 ierr = PetscOptionsReal("-tao_ntl_eta3", "good steplength; increase radius", "", tl->eta3, &tl->eta3,NULL);CHKERRQ(ierr); 655 ierr = PetscOptionsReal("-tao_ntl_eta4", "excellent steplength; greatly increase radius", "", tl->eta4, &tl->eta4,NULL);CHKERRQ(ierr); 656 ierr = PetscOptionsReal("-tao_ntl_alpha1", "", "", tl->alpha1, &tl->alpha1,NULL);CHKERRQ(ierr); 657 ierr = PetscOptionsReal("-tao_ntl_alpha2", "", "", tl->alpha2, &tl->alpha2,NULL);CHKERRQ(ierr); 658 ierr = PetscOptionsReal("-tao_ntl_alpha3", "", "", tl->alpha3, &tl->alpha3,NULL);CHKERRQ(ierr); 659 ierr = PetscOptionsReal("-tao_ntl_alpha4", "", "", tl->alpha4, &tl->alpha4,NULL);CHKERRQ(ierr); 660 ierr = PetscOptionsReal("-tao_ntl_alpha5", "", "", tl->alpha5, &tl->alpha5,NULL);CHKERRQ(ierr); 661 ierr = PetscOptionsReal("-tao_ntl_nu1", "poor steplength; reduce radius", "", tl->nu1, &tl->nu1,NULL);CHKERRQ(ierr); 662 ierr = PetscOptionsReal("-tao_ntl_nu2", "reasonable steplength; leave radius alone", "", tl->nu2, &tl->nu2,NULL);CHKERRQ(ierr); 663 ierr = PetscOptionsReal("-tao_ntl_nu3", "good steplength; increase radius", "", tl->nu3, &tl->nu3,NULL);CHKERRQ(ierr); 664 ierr = PetscOptionsReal("-tao_ntl_nu4", "excellent steplength; greatly increase radius", "", tl->nu4, &tl->nu4,NULL);CHKERRQ(ierr); 665 ierr = PetscOptionsReal("-tao_ntl_omega1", "", "", tl->omega1, &tl->omega1,NULL);CHKERRQ(ierr); 666 ierr = PetscOptionsReal("-tao_ntl_omega2", "", "", tl->omega2, &tl->omega2,NULL);CHKERRQ(ierr); 667 ierr = PetscOptionsReal("-tao_ntl_omega3", "", "", tl->omega3, &tl->omega3,NULL);CHKERRQ(ierr); 668 ierr = PetscOptionsReal("-tao_ntl_omega4", "", "", tl->omega4, &tl->omega4,NULL);CHKERRQ(ierr); 669 ierr = PetscOptionsReal("-tao_ntl_omega5", "", "", tl->omega5, &tl->omega5,NULL);CHKERRQ(ierr); 670 ierr = PetscOptionsReal("-tao_ntl_mu1_i", "", "", tl->mu1_i, &tl->mu1_i,NULL);CHKERRQ(ierr); 671 ierr = PetscOptionsReal("-tao_ntl_mu2_i", "", "", tl->mu2_i, &tl->mu2_i,NULL);CHKERRQ(ierr); 672 ierr = PetscOptionsReal("-tao_ntl_gamma1_i", "", "", tl->gamma1_i, &tl->gamma1_i,NULL);CHKERRQ(ierr); 673 ierr = PetscOptionsReal("-tao_ntl_gamma2_i", "", "", tl->gamma2_i, &tl->gamma2_i,NULL);CHKERRQ(ierr); 674 ierr = PetscOptionsReal("-tao_ntl_gamma3_i", "", "", tl->gamma3_i, &tl->gamma3_i,NULL);CHKERRQ(ierr); 675 ierr = PetscOptionsReal("-tao_ntl_gamma4_i", "", "", tl->gamma4_i, &tl->gamma4_i,NULL);CHKERRQ(ierr); 676 ierr = PetscOptionsReal("-tao_ntl_theta_i", "", "", tl->theta_i, &tl->theta_i,NULL);CHKERRQ(ierr); 677 ierr = PetscOptionsReal("-tao_ntl_mu1", "", "", tl->mu1, &tl->mu1,NULL);CHKERRQ(ierr); 678 ierr = PetscOptionsReal("-tao_ntl_mu2", "", "", tl->mu2, &tl->mu2,NULL);CHKERRQ(ierr); 679 ierr = PetscOptionsReal("-tao_ntl_gamma1", "", "", tl->gamma1, &tl->gamma1,NULL);CHKERRQ(ierr); 680 ierr = PetscOptionsReal("-tao_ntl_gamma2", "", "", tl->gamma2, &tl->gamma2,NULL);CHKERRQ(ierr); 681 ierr = PetscOptionsReal("-tao_ntl_gamma3", "", "", tl->gamma3, &tl->gamma3,NULL);CHKERRQ(ierr); 682 ierr = PetscOptionsReal("-tao_ntl_gamma4", "", "", tl->gamma4, &tl->gamma4,NULL);CHKERRQ(ierr); 683 ierr = PetscOptionsReal("-tao_ntl_theta", "", "", tl->theta, &tl->theta,NULL);CHKERRQ(ierr); 684 ierr = PetscOptionsReal("-tao_ntl_min_radius", "lower bound on initial radius", "", tl->min_radius, &tl->min_radius,NULL);CHKERRQ(ierr); 685 ierr = PetscOptionsReal("-tao_ntl_max_radius", "upper bound on radius", "", tl->max_radius, &tl->max_radius,NULL);CHKERRQ(ierr); 686 ierr = PetscOptionsReal("-tao_ntl_epsilon", "tolerance used when computing actual and predicted reduction", "", tl->epsilon, &tl->epsilon,NULL);CHKERRQ(ierr); 687 ierr = PetscOptionsTail();CHKERRQ(ierr); 688 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 689 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 690 PetscFunctionReturn(0); 691 } 692 693 /*------------------------------------------------------------*/ 694 static PetscErrorCode TaoView_NTL(Tao tao, PetscViewer viewer) 695 { 696 TAO_NTL *tl = (TAO_NTL *)tao->data; 697 PetscBool isascii; 698 PetscErrorCode ierr; 699 700 PetscFunctionBegin; 701 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 702 if (isascii) { 703 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 704 ierr = PetscViewerASCIIPrintf(viewer, "Trust-region steps: %D\n", tl->ntrust);CHKERRQ(ierr); 705 ierr = PetscViewerASCIIPrintf(viewer, "Newton search steps: %D\n", tl->newt);CHKERRQ(ierr); 706 ierr = PetscViewerASCIIPrintf(viewer, "BFGS search steps: %D\n", tl->bfgs);CHKERRQ(ierr); 707 ierr = PetscViewerASCIIPrintf(viewer, "Gradient search steps: %D\n", tl->grad);CHKERRQ(ierr); 708 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 709 } 710 PetscFunctionReturn(0); 711 } 712 713 /* ---------------------------------------------------------- */ 714 /*MC 715 TAONTR - Newton's method with trust region and linesearch 716 for unconstrained minimization. 717 At each iteration, the Newton trust region method solves the system for d 718 and performs a line search in the d direction: 719 720 min_d .5 dT Hk d + gkT d, s.t. ||d|| < Delta_k 721 722 Options Database Keys: 723 + -tao_ntl_pc_type - "none","ahess","bfgs","petsc" 724 . -tao_ntl_bfgs_scale_type - type of scaling with bfgs pc, "ahess" or "bfgs" 725 . -tao_ntl_init_type - "constant","direction","interpolation" 726 . -tao_ntl_update_type - "reduction","interpolation" 727 . -tao_ntl_min_radius - lower bound on trust region radius 728 . -tao_ntl_max_radius - upper bound on trust region radius 729 . -tao_ntl_epsilon - tolerance for accepting actual / predicted reduction 730 . -tao_ntl_mu1_i - mu1 interpolation init factor 731 . -tao_ntl_mu2_i - mu2 interpolation init factor 732 . -tao_ntl_gamma1_i - gamma1 interpolation init factor 733 . -tao_ntl_gamma2_i - gamma2 interpolation init factor 734 . -tao_ntl_gamma3_i - gamma3 interpolation init factor 735 . -tao_ntl_gamma4_i - gamma4 interpolation init factor 736 . -tao_ntl_theta_i - thetha1 interpolation init factor 737 . -tao_ntl_eta1 - eta1 reduction update factor 738 . -tao_ntl_eta2 - eta2 reduction update factor 739 . -tao_ntl_eta3 - eta3 reduction update factor 740 . -tao_ntl_eta4 - eta4 reduction update factor 741 . -tao_ntl_alpha1 - alpha1 reduction update factor 742 . -tao_ntl_alpha2 - alpha2 reduction update factor 743 . -tao_ntl_alpha3 - alpha3 reduction update factor 744 . -tao_ntl_alpha4 - alpha4 reduction update factor 745 . -tao_ntl_alpha4 - alpha4 reduction update factor 746 . -tao_ntl_mu1 - mu1 interpolation update 747 . -tao_ntl_mu2 - mu2 interpolation update 748 . -tao_ntl_gamma1 - gamma1 interpolcation update 749 . -tao_ntl_gamma2 - gamma2 interpolcation update 750 . -tao_ntl_gamma3 - gamma3 interpolcation update 751 . -tao_ntl_gamma4 - gamma4 interpolation update 752 - -tao_ntl_theta - theta1 interpolation update 753 754 Level: beginner 755 M*/ 756 757 PETSC_EXTERN PetscErrorCode TaoCreate_NTL(Tao tao) 758 { 759 TAO_NTL *tl; 760 PetscErrorCode ierr; 761 const char *morethuente_type = TAOLINESEARCHMT; 762 763 PetscFunctionBegin; 764 ierr = PetscNewLog(tao,&tl);CHKERRQ(ierr); 765 tao->ops->setup = TaoSetUp_NTL; 766 tao->ops->solve = TaoSolve_NTL; 767 tao->ops->view = TaoView_NTL; 768 tao->ops->setfromoptions = TaoSetFromOptions_NTL; 769 tao->ops->destroy = TaoDestroy_NTL; 770 771 /* Override default settings (unless already changed) */ 772 if (!tao->max_it_changed) tao->max_it = 50; 773 if (!tao->trust0_changed) tao->trust0 = 100.0; 774 775 tao->data = (void*)tl; 776 777 /* Default values for trust-region radius update based on steplength */ 778 tl->nu1 = 0.25; 779 tl->nu2 = 0.50; 780 tl->nu3 = 1.00; 781 tl->nu4 = 1.25; 782 783 tl->omega1 = 0.25; 784 tl->omega2 = 0.50; 785 tl->omega3 = 1.00; 786 tl->omega4 = 2.00; 787 tl->omega5 = 4.00; 788 789 /* Default values for trust-region radius update based on reduction */ 790 tl->eta1 = 1.0e-4; 791 tl->eta2 = 0.25; 792 tl->eta3 = 0.50; 793 tl->eta4 = 0.90; 794 795 tl->alpha1 = 0.25; 796 tl->alpha2 = 0.50; 797 tl->alpha3 = 1.00; 798 tl->alpha4 = 2.00; 799 tl->alpha5 = 4.00; 800 801 /* Default values for trust-region radius update based on interpolation */ 802 tl->mu1 = 0.10; 803 tl->mu2 = 0.50; 804 805 tl->gamma1 = 0.25; 806 tl->gamma2 = 0.50; 807 tl->gamma3 = 2.00; 808 tl->gamma4 = 4.00; 809 810 tl->theta = 0.05; 811 812 /* Default values for trust region initialization based on interpolation */ 813 tl->mu1_i = 0.35; 814 tl->mu2_i = 0.50; 815 816 tl->gamma1_i = 0.0625; 817 tl->gamma2_i = 0.5; 818 tl->gamma3_i = 2.0; 819 tl->gamma4_i = 5.0; 820 821 tl->theta_i = 0.25; 822 823 /* Remaining parameters */ 824 tl->min_radius = 1.0e-10; 825 tl->max_radius = 1.0e10; 826 tl->epsilon = 1.0e-6; 827 828 tl->init_type = NTL_INIT_INTERPOLATION; 829 tl->update_type = NTL_UPDATE_REDUCTION; 830 831 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 832 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 833 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 834 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 835 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 836 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 837 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 838 ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 839 ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842