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 = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 22 PetscErrorCode TaoBNKInitialize(Tao tao) 23 { 24 PetscErrorCode ierr; 25 TAO_BNK *bnk = (TAO_BNK *)tao->data; 26 KSPType ksp_type; 27 PC pc; 28 29 PetscReal fmin, ftrial, prered, actred, kappa, sigma; 30 PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 31 PetscReal delta, step = 1.0; 32 33 PetscInt n,N,needH = 1; 34 35 PetscInt i_max = 5; 36 PetscInt j_max = 1; 37 PetscInt i, j; 38 39 PetscFunctionBegin; 40 /* Number of times ksp stopped because of these reasons */ 41 bnk->ksp_atol = 0; 42 bnk->ksp_rtol = 0; 43 bnk->ksp_dtol = 0; 44 bnk->ksp_ctol = 0; 45 bnk->ksp_negc = 0; 46 bnk->ksp_iter = 0; 47 bnk->ksp_othr = 0; 48 49 /* Initialize trust-region radius when using nash, stcg, or gltr 50 Command automatically ignored for other methods 51 Will be reset during the first iteration 52 */ 53 ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 54 ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 55 ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 56 ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 57 58 ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 59 60 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 61 if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative"); 62 tao->trust = tao->trust0; 63 tao->trust = PetscMax(tao->trust, bnk->min_radius); 64 tao->trust = PetscMin(tao->trust, bnk->max_radius); 65 } 66 67 /* Get vectors we will need */ 68 if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) { 69 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 70 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 71 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr); 72 ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr); 73 } 74 75 /* create vectors for the limited memory preconditioner */ 76 if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) { 77 if (!bnk->Diag) { 78 ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr); 79 } 80 } 81 82 /* Modify the preconditioner to use the bfgs approximation */ 83 ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 84 switch(bnk->pc_type) { 85 case BNK_PC_NONE: 86 ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 87 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 88 break; 89 90 case BNK_PC_AHESS: 91 ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 92 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 93 ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 94 break; 95 96 case BNK_PC_BFGS: 97 ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 98 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 99 ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 100 ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr); 101 ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 102 break; 103 104 default: 105 /* Use the pc method set by pc_type */ 106 break; 107 } 108 109 /* Initialize trust-region radius. The initialization is only performed 110 when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 111 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 112 switch(bnk->init_type) { 113 case BNK_INIT_CONSTANT: 114 /* Use the initial radius specified */ 115 break; 116 117 case BNK_INIT_INTERPOLATION: 118 /* Use the initial radius specified */ 119 max_radius = 0.0; 120 121 for (j = 0; j < j_max; ++j) { 122 fmin = bnk->f; 123 sigma = 0.0; 124 125 if (needH) { 126 ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 127 needH = 0; 128 } 129 130 for (i = 0; i < i_max; ++i) { 131 ierr = VecCopy(tao->solution,bnk->W);CHKERRQ(ierr); 132 ierr = VecAXPY(bnk->W,-tao->trust/bnk->gnorm,tao->gradient);CHKERRQ(ierr); 133 ierr = TaoComputeObjective(tao, bnk->W, &ftrial);CHKERRQ(ierr); 134 if (PetscIsInfOrNanReal(ftrial)) { 135 tau = bnk->gamma1_i; 136 } else { 137 if (ftrial < fmin) { 138 fmin = ftrial; 139 sigma = -tao->trust / bnk->gnorm; 140 } 141 142 ierr = MatMult(tao->hessian, tao->gradient, bnk->D);CHKERRQ(ierr); 143 ierr = VecDot(tao->gradient, bnk->D, &prered);CHKERRQ(ierr); 144 145 prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 146 actred = bnk->f - ftrial; 147 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 148 kappa = 1.0; 149 } else { 150 kappa = actred / prered; 151 } 152 153 tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 154 tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 155 tau_min = PetscMin(tau_1, tau_2); 156 tau_max = PetscMax(tau_1, tau_2); 157 158 if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) { 159 /* Great agreement */ 160 max_radius = PetscMax(max_radius, tao->trust); 161 162 if (tau_max < 1.0) { 163 tau = bnk->gamma3_i; 164 } else if (tau_max > bnk->gamma4_i) { 165 tau = bnk->gamma4_i; 166 } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) { 167 tau = tau_1; 168 } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) { 169 tau = tau_2; 170 } else { 171 tau = tau_max; 172 } 173 } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) { 174 /* Good agreement */ 175 max_radius = PetscMax(max_radius, tao->trust); 176 177 if (tau_max < bnk->gamma2_i) { 178 tau = bnk->gamma2_i; 179 } else if (tau_max > bnk->gamma3_i) { 180 tau = bnk->gamma3_i; 181 } else { 182 tau = tau_max; 183 } 184 } else { 185 /* Not good agreement */ 186 if (tau_min > 1.0) { 187 tau = bnk->gamma2_i; 188 } else if (tau_max < bnk->gamma1_i) { 189 tau = bnk->gamma1_i; 190 } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 191 tau = bnk->gamma1_i; 192 } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 193 tau = tau_1; 194 } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 195 tau = tau_2; 196 } else { 197 tau = tau_max; 198 } 199 } 200 } 201 tao->trust = tau * tao->trust; 202 } 203 204 if (fmin < bnk->f) { 205 bnk->f = fmin; 206 ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 207 ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr); 208 209 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 210 if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 211 needH = 1; 212 213 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 214 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 215 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 216 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 217 } 218 } 219 tao->trust = PetscMax(tao->trust, max_radius); 220 221 /* Modify the radius if it is too large or small */ 222 tao->trust = PetscMax(tao->trust, bnk->min_radius); 223 tao->trust = PetscMin(tao->trust, bnk->max_radius); 224 break; 225 226 default: 227 /* Norm of the first direction will initialize radius */ 228 tao->trust = 0.0; 229 break; 230 } 231 } 232 233 /* Set initial scaling for the BFGS preconditioner 234 This step is done after computing the initial trust-region radius 235 since the function value may have decreased */ 236 if (BNK_PC_BFGS == bnk->pc_type) { 237 if (bnk->f != 0.0) { 238 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 239 } else { 240 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 241 } 242 ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 243 } 244 245 /* Set counter for gradient/reset steps*/ 246 bnk->newt = 0; 247 bnk->bfgs = 0; 248 bnk->sgrad = 0; 249 bnk->grad = 0; 250 PetscFunctionReturn(0); 251 } 252 253 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscInt *stepType) 254 { 255 PetscErrorCode ierr; 256 TAO_BNK *bnk = (TAO_BNK *)tao->data; 257 KSPConvergedReason ksp_reason; 258 259 PetscReal gdx, delta; 260 PetscReal norm_d = 0.0, e_min; 261 262 PetscInt bfgsUpdates = 0; 263 PetscInt kspits; 264 PetscInt needH = 1; 265 266 PetscFunctionBegin; 267 /* Compute the Hessian */ 268 if (needH) { 269 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 270 } 271 272 if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_AHESS == bnk->bfgs_scale_type)) { 273 /* Obtain diagonal for the bfgs preconditioner */ 274 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 275 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 276 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 277 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 278 } 279 280 /* Shift the Hessian matrix */ 281 bnk->pert = bnk->sval; 282 if (bnk->pert > 0) { 283 ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr); 284 if (tao->hessian != tao->hessian_pre) { 285 ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr); 286 } 287 } 288 289 if (BNK_PC_BFGS == bnk->pc_type) { 290 if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 291 /* Obtain diagonal for the bfgs preconditioner */ 292 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 293 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 294 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 295 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 296 } 297 /* Update the limited memory preconditioner and get existing # of updates */ 298 ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 299 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 300 } 301 302 /* Solve the Newton system of equations */ 303 ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 304 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 305 ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 306 ierr = KSPSolve(tao->ksp, tao->gradient, bnk->D);CHKERRQ(ierr); 307 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 308 tao->ksp_its+=kspits; 309 tao->ksp_tot_its+=kspits; 310 ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 311 312 if (0.0 == tao->trust) { 313 /* Radius was uninitialized; use the norm of the direction */ 314 if (norm_d > 0.0) { 315 tao->trust = norm_d; 316 317 /* Modify the radius if it is too large or small */ 318 tao->trust = PetscMax(tao->trust, bnk->min_radius); 319 tao->trust = PetscMin(tao->trust, bnk->max_radius); 320 } else { 321 /* The direction was bad; set radius to default value and re-solve 322 the trust-region subproblem to get a direction */ 323 tao->trust = tao->trust0; 324 325 /* Modify the radius if it is too large or small */ 326 tao->trust = PetscMax(tao->trust, bnk->min_radius); 327 tao->trust = PetscMin(tao->trust, bnk->max_radius); 328 329 ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 330 ierr = KSPSolve(tao->ksp, tao->gradient, bnk->D);CHKERRQ(ierr); 331 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 332 tao->ksp_its+=kspits; 333 tao->ksp_tot_its+=kspits; 334 ierr = KSPCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 335 336 if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 337 } 338 } 339 } else { 340 ierr = KSPSolve(tao->ksp, tao->gradient, bnk->D);CHKERRQ(ierr); 341 ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 342 tao->ksp_its += kspits; 343 tao->ksp_tot_its+=kspits; 344 } 345 ierr = VecScale(bnk->D, -1.0);CHKERRQ(ierr); 346 ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 347 if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (BNK_PC_BFGS == bnk->pc_type) && (bfgsUpdates > 1)) { 348 /* Preconditioner is numerically indefinite; reset the 349 approximate if using BFGS preconditioning. */ 350 351 if (bnk->f != 0.0) { 352 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 353 } else { 354 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 355 } 356 ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 357 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 358 ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 359 bfgsUpdates = 1; 360 } 361 362 if (KSP_CONVERGED_ATOL == ksp_reason) { 363 ++bnk->ksp_atol; 364 } else if (KSP_CONVERGED_RTOL == ksp_reason) { 365 ++bnk->ksp_rtol; 366 } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 367 ++bnk->ksp_ctol; 368 } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 369 ++bnk->ksp_negc; 370 } else if (KSP_DIVERGED_DTOL == ksp_reason) { 371 ++bnk->ksp_dtol; 372 } else if (KSP_DIVERGED_ITS == ksp_reason) { 373 ++bnk->ksp_iter; 374 } else { 375 ++bnk->ksp_othr; 376 } 377 378 /* Check for success (descent direction) */ 379 ierr = VecDot(bnk->D, tao->gradient, &gdx);CHKERRQ(ierr); 380 if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 381 /* Newton step is not descent or direction produced Inf or NaN 382 Update the perturbation for next time */ 383 if (bnk->pert <= 0.0) { 384 /* Initialize the perturbation */ 385 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 386 if (bnk->is_gltr) { 387 ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 388 bnk->pert = PetscMax(bnk->pert, -e_min); 389 } 390 } else { 391 /* Increase the perturbation */ 392 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 393 } 394 395 if (BNK_PC_BFGS != bnk->pc_type) { 396 /* We don't have the bfgs matrix around and updated 397 Must use gradient direction in this case */ 398 ierr = VecCopy(tao->gradient, bnk->D);CHKERRQ(ierr); 399 ierr = VecScale(bnk->D, -1.0);CHKERRQ(ierr); 400 ++bnk->grad; 401 *stepType = BNK_GRADIENT; 402 } else { 403 /* Attempt to use the BFGS direction */ 404 ierr = MatLMVMSolve(bnk->M, tao->gradient, bnk->D);CHKERRQ(ierr); 405 ierr = VecScale(bnk->D, -1.0);CHKERRQ(ierr); 406 407 /* Check for success (descent direction) */ 408 ierr = VecDot(tao->gradient, bnk->D, &gdx);CHKERRQ(ierr); 409 if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 410 /* BFGS direction is not descent or direction produced not a number 411 We can assert bfgsUpdates > 1 in this case because 412 the first solve produces the scaled gradient direction, 413 which is guaranteed to be descent */ 414 415 /* Use steepest descent direction (scaled) */ 416 417 if (bnk->f != 0.0) { 418 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 419 } else { 420 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 421 } 422 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 423 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 424 ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 425 ierr = MatLMVMSolve(bnk->M, tao->gradient, bnk->D);CHKERRQ(ierr); 426 ierr = VecScale(bnk->D, -1.0);CHKERRQ(ierr); 427 428 bfgsUpdates = 1; 429 ++bnk->sgrad; 430 *stepType = BNK_SCALED_GRADIENT; 431 } else { 432 if (1 == bfgsUpdates) { 433 /* The first BFGS direction is always the scaled gradient */ 434 ++bnk->sgrad; 435 *stepType = BNK_SCALED_GRADIENT; 436 } else { 437 ++bnk->bfgs; 438 *stepType = BNK_BFGS; 439 } 440 } 441 } 442 } else { 443 /* Computed Newton step is descent */ 444 switch (ksp_reason) { 445 case KSP_DIVERGED_NANORINF: 446 case KSP_DIVERGED_BREAKDOWN: 447 case KSP_DIVERGED_INDEFINITE_MAT: 448 case KSP_DIVERGED_INDEFINITE_PC: 449 case KSP_CONVERGED_CG_NEG_CURVE: 450 /* Matrix or preconditioner is indefinite; increase perturbation */ 451 if (bnk->pert <= 0.0) { 452 /* Initialize the perturbation */ 453 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 454 if (bnk->is_gltr) { 455 ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 456 bnk->pert = PetscMax(bnk->pert, -e_min); 457 } 458 } else { 459 /* Increase the perturbation */ 460 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 461 } 462 break; 463 464 default: 465 /* Newton step computation is good; decrease perturbation */ 466 bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 467 if (bnk->pert < bnk->pmin) { 468 bnk->pert = 0.0; 469 } 470 break; 471 } 472 473 ++bnk->newt; 474 stepType = BNK_NEWTON; 475 } 476 PetscFunctionReturn(0); 477 } 478 479 /* ---------------------------------------------------------- */ 480 static PetscErrorCode TaoSetUp_BNK(Tao tao) 481 { 482 TAO_BNK *bnk = (TAO_BNK *)tao->data; 483 PetscErrorCode ierr; 484 485 PetscFunctionBegin; 486 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 487 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 488 if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 489 if (!bnk->D) {ierr = VecDuplicate(tao->solution,&bnk->D);CHKERRQ(ierr);} 490 if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 491 if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 492 bnk->Diag = 0; 493 bnk->M = 0; 494 PetscFunctionReturn(0); 495 } 496 497 /*------------------------------------------------------------*/ 498 static PetscErrorCode TaoDestroy_BNK(Tao tao) 499 { 500 TAO_BNK *bnk = (TAO_BNK *)tao->data; 501 PetscErrorCode ierr; 502 503 PetscFunctionBegin; 504 if (tao->setupcalled) { 505 ierr = VecDestroy(&bnk->D);CHKERRQ(ierr); 506 ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 507 ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 508 ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 509 } 510 ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 511 ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 512 ierr = PetscFree(tao->data);CHKERRQ(ierr); 513 PetscFunctionReturn(0); 514 } 515 516 /*------------------------------------------------------------*/ 517 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 518 { 519 TAO_BNK *bnk = (TAO_BNK *)tao->data; 520 PetscErrorCode ierr; 521 522 PetscFunctionBegin; 523 ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 524 ierr = PetscOptionsEList("-tao_BNK_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr); 525 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); 526 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); 527 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); 528 ierr = PetscOptionsReal("-tao_BNK_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 529 ierr = PetscOptionsReal("-tao_BNK_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 530 ierr = PetscOptionsReal("-tao_BNK_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 531 ierr = PetscOptionsReal("-tao_BNK_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 532 ierr = PetscOptionsReal("-tao_BNK_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 533 ierr = PetscOptionsReal("-tao_BNK_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 534 ierr = PetscOptionsReal("-tao_BNK_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 535 ierr = PetscOptionsReal("-tao_BNK_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 536 ierr = PetscOptionsReal("-tao_BNK_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 537 ierr = PetscOptionsReal("-tao_BNK_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 538 ierr = PetscOptionsReal("-tao_BNK_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 539 ierr = PetscOptionsReal("-tao_BNK_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 540 ierr = PetscOptionsReal("-tao_BNK_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 541 ierr = PetscOptionsReal("-tao_BNK_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 542 ierr = PetscOptionsReal("-tao_BNK_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 543 ierr = PetscOptionsReal("-tao_BNK_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 544 ierr = PetscOptionsReal("-tao_BNK_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 545 ierr = PetscOptionsReal("-tao_BNK_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 546 ierr = PetscOptionsReal("-tao_BNK_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 547 ierr = PetscOptionsReal("-tao_BNK_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 548 ierr = PetscOptionsReal("-tao_BNK_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 549 ierr = PetscOptionsReal("-tao_BNK_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 550 ierr = PetscOptionsReal("-tao_BNK_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 551 ierr = PetscOptionsReal("-tao_BNK_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 552 ierr = PetscOptionsReal("-tao_BNK_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 553 ierr = PetscOptionsReal("-tao_BNK_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 554 ierr = PetscOptionsReal("-tao_BNK_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 555 ierr = PetscOptionsReal("-tao_BNK_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 556 ierr = PetscOptionsReal("-tao_BNK_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 557 ierr = PetscOptionsReal("-tao_BNK_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 558 ierr = PetscOptionsReal("-tao_BNK_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 559 ierr = PetscOptionsReal("-tao_BNK_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 560 ierr = PetscOptionsReal("-tao_BNK_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 561 ierr = PetscOptionsReal("-tao_BNK_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 562 ierr = PetscOptionsReal("-tao_BNK_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 563 ierr = PetscOptionsReal("-tao_BNK_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 564 ierr = PetscOptionsReal("-tao_BNK_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 565 ierr = PetscOptionsReal("-tao_BNK_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 566 ierr = PetscOptionsReal("-tao_BNK_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 567 ierr = PetscOptionsReal("-tao_BNK_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 568 ierr = PetscOptionsReal("-tao_BNK_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 569 ierr = PetscOptionsReal("-tao_BNK_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 570 ierr = PetscOptionsReal("-tao_BNK_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 571 ierr = PetscOptionsReal("-tao_BNK_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 572 ierr = PetscOptionsReal("-tao_BNK_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 573 ierr = PetscOptionsTail();CHKERRQ(ierr); 574 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 575 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 580 /*------------------------------------------------------------*/ 581 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 582 { 583 TAO_BNK *bnk = (TAO_BNK *)tao->data; 584 PetscInt nrejects; 585 PetscBool isascii; 586 PetscErrorCode ierr; 587 588 PetscFunctionBegin; 589 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 590 if (isascii) { 591 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 592 if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 593 ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 594 ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 595 } 596 ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 597 ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 598 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 599 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 600 ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 601 ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 602 ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 603 ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 604 ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 605 ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 606 ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 607 ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 608 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 609 } 610 PetscFunctionReturn(0); 611 } 612 613 /* ---------------------------------------------------------- */ 614 /*MC 615 TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 616 At each iteration, the BNK method solves the symmetric 617 system of equations to obtain the step diretion dk: 618 Hk dk = -gk 619 at which point the step can be globalized either through trust-region 620 methods, or a line search, or a heuristic mixture of both. 621 622 Options Database Keys: 623 + -tao_BNK_pc_type - "none","ahess","bfgs","petsc" 624 . -tao_BNK_bfgs_scale_type - "ahess","phess","bfgs" 625 . -tao_BNK_init_type - "constant","direction","interpolation" 626 . -tao_BNK_update_type - "step","direction","interpolation" 627 . -tao_BNK_sval - perturbation starting value 628 . -tao_BNK_imin - minimum initial perturbation 629 . -tao_BNK_imax - maximum initial perturbation 630 . -tao_BNK_pmin - minimum perturbation 631 . -tao_BNK_pmax - maximum perturbation 632 . -tao_BNK_pgfac - growth factor 633 . -tao_BNK_psfac - shrink factor 634 . -tao_BNK_imfac - initial merit factor 635 . -tao_BNK_pmgfac - merit growth factor 636 . -tao_BNK_pmsfac - merit shrink factor 637 . -tao_BNK_eta1 - poor steplength; reduce radius 638 . -tao_BNK_eta2 - reasonable steplength; leave radius 639 . -tao_BNK_eta3 - good steplength; increase readius 640 . -tao_BNK_eta4 - excellent steplength; greatly increase radius 641 . -tao_BNK_alpha1 - alpha1 reduction 642 . -tao_BNK_alpha2 - alpha2 reduction 643 . -tao_BNK_alpha3 - alpha3 reduction 644 . -tao_BNK_alpha4 - alpha4 reduction 645 . -tao_BNK_alpha - alpha5 reduction 646 . -tao_BNK_mu1 - mu1 interpolation update 647 . -tao_BNK_mu2 - mu2 interpolation update 648 . -tao_BNK_gamma1 - gamma1 interpolation update 649 . -tao_BNK_gamma2 - gamma2 interpolation update 650 . -tao_BNK_gamma3 - gamma3 interpolation update 651 . -tao_BNK_gamma4 - gamma4 interpolation update 652 . -tao_BNK_theta - theta interpolation update 653 . -tao_BNK_omega1 - omega1 step update 654 . -tao_BNK_omega2 - omega2 step update 655 . -tao_BNK_omega3 - omega3 step update 656 . -tao_BNK_omega4 - omega4 step update 657 . -tao_BNK_omega5 - omega5 step update 658 . -tao_BNK_mu1_i - mu1 interpolation init factor 659 . -tao_BNK_mu2_i - mu2 interpolation init factor 660 . -tao_BNK_gamma1_i - gamma1 interpolation init factor 661 . -tao_BNK_gamma2_i - gamma2 interpolation init factor 662 . -tao_BNK_gamma3_i - gamma3 interpolation init factor 663 . -tao_BNK_gamma4_i - gamma4 interpolation init factor 664 - -tao_BNK_theta_i - theta interpolation init factor 665 666 Level: beginner 667 M*/ 668 669 PetscErrorCode TaoCreate_BNK(Tao tao) 670 { 671 TAO_BNK *bnk; 672 const char *morethuente_type = TAOLINESEARCHMT; 673 PetscErrorCode ierr; 674 675 PetscFunctionBegin; 676 ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 677 678 tao->ops->setup = TaoSetUp_BNK; 679 tao->ops->view = TaoView_BNK; 680 tao->ops->setfromoptions = TaoSetFromOptions_BNK; 681 tao->ops->destroy = TaoDestroy_BNK; 682 683 /* Override default settings (unless already changed) */ 684 if (!tao->max_it_changed) tao->max_it = 50; 685 if (!tao->trust0_changed) tao->trust0 = 100.0; 686 687 tao->data = (void*)bnk; 688 689 bnk->sval = 0.0; 690 bnk->imin = 1.0e-4; 691 bnk->imax = 1.0e+2; 692 bnk->imfac = 1.0e-1; 693 694 bnk->pmin = 1.0e-12; 695 bnk->pmax = 1.0e+2; 696 bnk->pgfac = 1.0e+1; 697 bnk->psfac = 4.0e-1; 698 bnk->pmgfac = 1.0e-1; 699 bnk->pmsfac = 1.0e-1; 700 701 /* Default values for trust-region radius update based on steplength */ 702 bnk->nu1 = 0.25; 703 bnk->nu2 = 0.50; 704 bnk->nu3 = 1.00; 705 bnk->nu4 = 1.25; 706 707 bnk->omega1 = 0.25; 708 bnk->omega2 = 0.50; 709 bnk->omega3 = 1.00; 710 bnk->omega4 = 2.00; 711 bnk->omega5 = 4.00; 712 713 /* Default values for trust-region radius update based on reduction */ 714 bnk->eta1 = 1.0e-4; 715 bnk->eta2 = 0.25; 716 bnk->eta3 = 0.50; 717 bnk->eta4 = 0.90; 718 719 bnk->alpha1 = 0.25; 720 bnk->alpha2 = 0.50; 721 bnk->alpha3 = 1.00; 722 bnk->alpha4 = 2.00; 723 bnk->alpha5 = 4.00; 724 725 /* Default values for trust-region radius update based on interpolation */ 726 bnk->mu1 = 0.10; 727 bnk->mu2 = 0.50; 728 729 bnk->gamma1 = 0.25; 730 bnk->gamma2 = 0.50; 731 bnk->gamma3 = 2.00; 732 bnk->gamma4 = 4.00; 733 734 bnk->theta = 0.05; 735 736 /* Default values for trust region initialization based on interpolation */ 737 bnk->mu1_i = 0.35; 738 bnk->mu2_i = 0.50; 739 740 bnk->gamma1_i = 0.0625; 741 bnk->gamma2_i = 0.5; 742 bnk->gamma3_i = 2.0; 743 bnk->gamma4_i = 5.0; 744 745 bnk->theta_i = 0.25; 746 747 /* Remaining parameters */ 748 bnk->min_radius = 1.0e-10; 749 bnk->max_radius = 1.0e10; 750 bnk->epsilon = 1.0e-6; 751 752 bnk->pc_type = BNK_PC_BFGS; 753 bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 754 bnk->init_type = BNK_INIT_INTERPOLATION; 755 bnk->update_type = BNK_UPDATE_STEP; 756 757 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 758 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 759 ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 760 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 761 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 762 763 /* Set linear solver to default for symmetric matrices */ 764 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 765 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 766 ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 767 ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 768 PetscFunctionReturn(0); 769 } 770