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