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 /*------------------------------------------------------------*/ 23 24 /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */ 25 26 PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType) 27 { 28 PetscErrorCode ierr; 29 TAO_BNK *bnk = (TAO_BNK *)tao->data; 30 PC pc; 31 32 PetscReal f_min, ftrial, prered, actred, kappa, sigma; 33 PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 34 PetscReal delta; 35 36 PetscInt n,N,needH = 1; 37 38 PetscInt i_max = 5; 39 PetscInt j_max = 1; 40 PetscInt i, j; 41 42 PetscFunctionBegin; 43 /* Project the current point onto the feasible set */ 44 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 45 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 46 47 /* Project the initial point onto the feasible region */ 48 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 49 50 /* Check convergence criteria */ 51 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 52 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 53 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 54 if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 55 56 /* Number of times ksp stopped because of these reasons */ 57 bnk->ksp_atol = 0; 58 bnk->ksp_rtol = 0; 59 bnk->ksp_dtol = 0; 60 bnk->ksp_ctol = 0; 61 bnk->ksp_negc = 0; 62 bnk->ksp_iter = 0; 63 bnk->ksp_othr = 0; 64 65 /* Initialize the Hessian perturbation */ 66 bnk->pert = bnk->sval; 67 68 /* Initialize trust-region radius when using nash, stcg, or gltr 69 Command automatically ignored for other methods 70 Will be reset during the first iteration 71 */ 72 ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr); 73 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 74 if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative"); 75 tao->trust = tao->trust0; 76 tao->trust = PetscMax(tao->trust, bnk->min_radius); 77 tao->trust = PetscMin(tao->trust, bnk->max_radius); 78 } 79 80 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 81 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr); 82 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 83 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 84 85 /* Get vectors we will need */ 86 if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) { 87 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 88 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 89 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr); 90 ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr); 91 } 92 93 /* create vectors for the limited memory preconditioner */ 94 if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) { 95 if (!bnk->Diag) {ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);} 96 } 97 ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr); 98 ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr); 99 100 /* Modify the preconditioner to use the bfgs approximation */ 101 ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 102 switch(bnk->pc_type) { 103 case BNK_PC_NONE: 104 ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 105 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 106 break; 107 108 case BNK_PC_AHESS: 109 ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 110 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 111 ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 112 break; 113 114 case BNK_PC_BFGS: 115 ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 116 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 117 ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 118 ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr); 119 ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 120 break; 121 122 default: 123 /* Use the pc method set by pc_type */ 124 break; 125 } 126 127 /* Initialize trust-region radius. The initialization is only performed 128 when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 129 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 130 switch(initType) { 131 case BNK_INIT_CONSTANT: 132 /* Use the initial radius specified */ 133 break; 134 135 case BNK_INIT_INTERPOLATION: 136 /* Use the initial radius specified */ 137 max_radius = 0.0; 138 139 for (j = 0; j < j_max; ++j) { 140 f_min = bnk->f; 141 sigma = 0.0; 142 143 if (needH) { 144 /* Compute the Hessian at the new step, and extract the inactive subsystem */ 145 ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 146 if (bnk->inactive_idx) { 147 ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 148 ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr); 149 ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr); 150 } else { 151 ierr = MatDestroy(&bnk->H_inactive); 152 ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive); 153 } 154 needH = 0; 155 } 156 157 for (i = 0; i < i_max; ++i) { 158 /* Take a steepest descent step and snap it to bounds */ 159 ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 160 ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr); 161 ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 162 /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */ 163 ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr); 164 ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr); 165 /* Compute the objective at the trial */ 166 ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr); 167 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 168 if (PetscIsInfOrNanReal(ftrial)) { 169 tau = bnk->gamma1_i; 170 } else { 171 if (ftrial < f_min) { 172 f_min = ftrial; 173 sigma = -tao->trust / bnk->gnorm; 174 } 175 /* Compute the predicted and actual reduction */ 176 ierr = MatMult(bnk->H_inactive, tao->gradient, bnk->W);CHKERRQ(ierr); 177 ierr = VecDot(tao->gradient, bnk->W, &prered);CHKERRQ(ierr); 178 prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm)); 179 actred = bnk->f - ftrial; 180 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 181 kappa = 1.0; 182 } else { 183 kappa = actred / prered; 184 } 185 186 tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred); 187 tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred); 188 tau_min = PetscMin(tau_1, tau_2); 189 tau_max = PetscMax(tau_1, tau_2); 190 191 if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) { 192 /* Great agreement */ 193 max_radius = PetscMax(max_radius, tao->trust); 194 195 if (tau_max < 1.0) { 196 tau = bnk->gamma3_i; 197 } else if (tau_max > bnk->gamma4_i) { 198 tau = bnk->gamma4_i; 199 } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) { 200 tau = tau_1; 201 } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) { 202 tau = tau_2; 203 } else { 204 tau = tau_max; 205 } 206 } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) { 207 /* Good agreement */ 208 max_radius = PetscMax(max_radius, tao->trust); 209 210 if (tau_max < bnk->gamma2_i) { 211 tau = bnk->gamma2_i; 212 } else if (tau_max > bnk->gamma3_i) { 213 tau = bnk->gamma3_i; 214 } else { 215 tau = tau_max; 216 } 217 } else { 218 /* Not good agreement */ 219 if (tau_min > 1.0) { 220 tau = bnk->gamma2_i; 221 } else if (tau_max < bnk->gamma1_i) { 222 tau = bnk->gamma1_i; 223 } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) { 224 tau = bnk->gamma1_i; 225 } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 226 tau = tau_1; 227 } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) { 228 tau = tau_2; 229 } else { 230 tau = tau_max; 231 } 232 } 233 } 234 tao->trust = tau * tao->trust; 235 } 236 237 if (f_min < bnk->f) { 238 bnk->f = f_min; 239 ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 240 ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 241 ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 242 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 243 244 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 245 if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 246 needH = 1; 247 248 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 249 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr); 250 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 251 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 252 } 253 } 254 tao->trust = PetscMax(tao->trust, max_radius); 255 256 /* Modify the radius if it is too large or small */ 257 tao->trust = PetscMax(tao->trust, bnk->min_radius); 258 tao->trust = PetscMin(tao->trust, bnk->max_radius); 259 break; 260 261 default: 262 /* Norm of the first direction will initialize radius */ 263 tao->trust = 0.0; 264 break; 265 } 266 } 267 268 /* Set initial scaling for the BFGS preconditioner 269 This step is done after computing the initial trust-region radius 270 since the function value may have decreased */ 271 if (BNK_PC_BFGS == bnk->pc_type) { 272 delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 273 ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 274 } 275 276 /* Set counter for gradient/reset steps*/ 277 bnk->newt = 0; 278 bnk->bfgs = 0; 279 bnk->sgrad = 0; 280 bnk->grad = 0; 281 PetscFunctionReturn(0); 282 } 283 284 /*------------------------------------------------------------*/ 285 286 /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */ 287 288 PetscErrorCode TaoBNKComputeHessian(Tao tao) 289 { 290 PetscErrorCode ierr; 291 TAO_BNK *bnk = (TAO_BNK *)tao->data; 292 293 PetscFunctionBegin; 294 /* Compute the Hessian */ 295 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 296 /* Add a correction to the BFGS preconditioner */ 297 if (BNK_PC_BFGS == bnk->pc_type) { 298 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 299 /* Update the BFGS diagonal scaling */ 300 if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type) { 301 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 302 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 303 ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr); 304 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 305 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 306 } 307 } 308 PetscFunctionReturn(0); 309 } 310 311 /*------------------------------------------------------------*/ 312 313 /* Routine for estimating the active set */ 314 315 PetscErrorCode TaoBNKEstimateActiveSet(Tao tao) 316 { 317 PetscErrorCode ierr; 318 TAO_BNK *bnk = (TAO_BNK *)tao->data; 319 320 PetscFunctionBegin; 321 switch (bnk->as_type) { 322 case BNK_AS_NONE: 323 ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr); 324 ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr); 325 ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr); 326 ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr); 327 break; 328 329 case BNK_AS_BERTSEKAS: 330 /* Compute the trial step vector with which we will estimate the active set at the next iteration */ 331 if (BNK_PC_BFGS == bnk->pc_type) { 332 /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */ 333 ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr); 334 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr); 335 } else { 336 /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */ 337 ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr); 338 ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr); 339 ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr); 340 ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr); 341 ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr); 342 } 343 ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr); 344 ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr); 345 346 default: 347 break; 348 } 349 PetscFunctionReturn(0); 350 } 351 352 /*------------------------------------------------------------*/ 353 354 /* Routine for bounding the step direction */ 355 356 PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step) 357 { 358 PetscErrorCode ierr; 359 TAO_BNK *bnk = (TAO_BNK *)tao->data; 360 361 PetscFunctionBegin; 362 if (bnk->active_idx) { 363 switch (bnk->as_type) { 364 case BNK_AS_NONE: 365 if (bnk->active_idx) { 366 ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 367 ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr); 368 ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr); 369 } 370 break; 371 372 case BNK_AS_BERTSEKAS: 373 ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr); 374 break; 375 376 default: 377 break; 378 } 379 } 380 PetscFunctionReturn(0); 381 } 382 383 /*------------------------------------------------------------*/ 384 385 /* Routine for computing the Newton step. 386 387 If the safeguard is enabled, the Newton step is verified to be a 388 descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled 389 gradient steps if/when necessary. 390 391 The function reports back on which type of step has ultimately been stored 392 under tao->stepdirection. 393 */ 394 395 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason) 396 { 397 PetscErrorCode ierr; 398 TAO_BNK *bnk = (TAO_BNK *)tao->data; 399 400 PetscReal delta; 401 PetscInt bfgsUpdates = 0; 402 PetscInt kspits; 403 404 PetscFunctionBegin; 405 /* Determine the active and inactive sets */ 406 ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr); 407 408 /* Prepare the reduced sub-matrices for the inactive set */ 409 if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); } 410 if (bnk->inactive_idx) { 411 ierr = MatDestroy(&bnk->H_inactive); 412 ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr); 413 if (tao->hessian == tao->hessian_pre) { 414 bnk->Hpre_inactive = bnk->H_inactive; 415 } else { 416 ierr = MatDestroy(&bnk->Hpre_inactive); 417 ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr); 418 } 419 } else { 420 ierr = MatDestroy(&bnk->H_inactive); 421 ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive); 422 if (tao->hessian == tao->hessian_pre) { 423 bnk->Hpre_inactive = bnk->H_inactive; 424 } else { 425 ierr = MatDestroy(&bnk->Hpre_inactive); 426 ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive); 427 } 428 } 429 430 /* Shift the reduced Hessian matrix */ 431 if ((shift) && (bnk->pert > 0)) { 432 ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr); 433 if (bnk->H_inactive != bnk->Hpre_inactive) { 434 ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr); 435 } 436 } 437 438 /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */ 439 if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) { 440 /* Obtain diagonal for the bfgs preconditioner */ 441 ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr); 442 if (bnk->inactive_idx) { 443 ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr); 444 } else { 445 bnk->Diag_red = bnk->Diag; 446 } 447 ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr); 448 if (bnk->inactive_idx) { 449 ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr); 450 } 451 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 452 ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr); 453 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 454 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 455 } 456 457 /* Solve the Newton system of equations */ 458 ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr); 459 ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 460 ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr); 461 ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr); 462 if (bnk->inactive_idx) { 463 ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 464 ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 465 } else { 466 bnk->G_inactive = bnk->unprojected_gradient; 467 bnk->X_inactive = tao->stepdirection; 468 } 469 if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) { 470 ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 471 ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr); 472 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 473 tao->ksp_its+=kspits; 474 tao->ksp_tot_its+=kspits; 475 ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 476 477 if (0.0 == tao->trust) { 478 /* Radius was uninitialized; use the norm of the direction */ 479 if (bnk->dnorm > 0.0) { 480 tao->trust = bnk->dnorm; 481 482 /* Modify the radius if it is too large or small */ 483 tao->trust = PetscMax(tao->trust, bnk->min_radius); 484 tao->trust = PetscMin(tao->trust, bnk->max_radius); 485 } else { 486 /* The direction was bad; set radius to default value and re-solve 487 the trust-region subproblem to get a direction */ 488 tao->trust = tao->trust0; 489 490 /* Modify the radius if it is too large or small */ 491 tao->trust = PetscMax(tao->trust, bnk->min_radius); 492 tao->trust = PetscMin(tao->trust, bnk->max_radius); 493 494 ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr); 495 ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr); 496 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 497 tao->ksp_its+=kspits; 498 tao->ksp_tot_its+=kspits; 499 ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr); 500 501 if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 502 } 503 } 504 } else { 505 ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr); 506 ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 507 tao->ksp_its += kspits; 508 tao->ksp_tot_its+=kspits; 509 } 510 /* Restore sub vectors back */ 511 if (bnk->inactive_idx) { 512 ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 513 ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 514 } 515 /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 516 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 517 ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 518 519 /* Record convergence reasons */ 520 ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr); 521 if (KSP_CONVERGED_ATOL == *ksp_reason) { 522 ++bnk->ksp_atol; 523 } else if (KSP_CONVERGED_RTOL == *ksp_reason) { 524 ++bnk->ksp_rtol; 525 } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) { 526 ++bnk->ksp_ctol; 527 } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) { 528 ++bnk->ksp_negc; 529 } else if (KSP_DIVERGED_DTOL == *ksp_reason) { 530 ++bnk->ksp_dtol; 531 } else if (KSP_DIVERGED_ITS == *ksp_reason) { 532 ++bnk->ksp_iter; 533 } else { 534 ++bnk->ksp_othr; 535 } 536 537 /* Make sure the BFGS preconditioner is healthy */ 538 if (bnk->pc_type == BNK_PC_BFGS) { 539 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 540 if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) { 541 /* Preconditioner is numerically indefinite; reset the approximation. */ 542 delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 543 ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr); 544 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 545 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 546 } 547 } 548 PetscFunctionReturn(0); 549 } 550 551 /*------------------------------------------------------------*/ 552 553 /* Routine for recomputing the predicted reduction for a given step vector */ 554 555 PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered) 556 { 557 PetscErrorCode ierr; 558 TAO_BNK *bnk = (TAO_BNK *)tao->data; 559 560 PetscFunctionBegin; 561 /* Extract subvectors associated with the inactive set */ 562 if (bnk->inactive_idx){ 563 ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 564 ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr); 565 ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 566 } else { 567 bnk->X_inactive = tao->stepdirection; 568 bnk->inactive_work = bnk->Xwork; 569 bnk->G_inactive = bnk->Gwork; 570 } 571 /* Recompute the predicted decrease based on the quadratic model */ 572 ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr); 573 ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr); 574 ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered); 575 /* Restore the sub vectors */ 576 if (bnk->inactive_idx){ 577 ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr); 578 ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr); 579 ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr); 580 } 581 PetscFunctionReturn(0); 582 } 583 584 /*------------------------------------------------------------*/ 585 586 /* Routine for ensuring that the Newton step is a descent direction. 587 588 The step direction falls back onto BFGS, scaled gradient and gradient steps 589 in the event that the Newton step fails the test. 590 */ 591 592 PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType) 593 { 594 PetscErrorCode ierr; 595 TAO_BNK *bnk = (TAO_BNK *)tao->data; 596 597 PetscReal gdx, delta, e_min; 598 PetscInt bfgsUpdates; 599 600 PetscFunctionBegin; 601 ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 602 if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 603 /* Newton step is not descent or direction produced Inf or NaN 604 Update the perturbation for next time */ 605 if (bnk->pert <= 0.0) { 606 /* Initialize the perturbation */ 607 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 608 if (bnk->is_gltr) { 609 ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 610 bnk->pert = PetscMax(bnk->pert, -e_min); 611 } 612 } else { 613 /* Increase the perturbation */ 614 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 615 } 616 617 if (BNK_PC_BFGS != bnk->pc_type) { 618 /* We don't have the bfgs matrix around and updated 619 Must use gradient direction in this case */ 620 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 621 *stepType = BNK_GRADIENT; 622 } else { 623 /* Attempt to use the BFGS direction */ 624 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 625 626 /* Check for success (descent direction) 627 NOTE: Negative gdx here means not a descent direction because 628 the fall-back step is missing a negative sign. */ 629 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 630 if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 631 /* BFGS direction is not descent or direction produced not a number 632 We can assert bfgsUpdates > 1 in this case because 633 the first solve produces the scaled gradient direction, 634 which is guaranteed to be descent */ 635 636 /* Use steepest descent direction (scaled) */ 637 delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 638 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 639 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 640 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 641 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 642 643 *stepType = BNK_SCALED_GRADIENT; 644 } else { 645 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 646 if (1 == bfgsUpdates) { 647 /* The first BFGS direction is always the scaled gradient */ 648 *stepType = BNK_SCALED_GRADIENT; 649 } else { 650 *stepType = BNK_BFGS; 651 } 652 } 653 } 654 /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 655 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 656 ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 657 } else { 658 /* Computed Newton step is descent */ 659 switch (ksp_reason) { 660 case KSP_DIVERGED_NANORINF: 661 case KSP_DIVERGED_BREAKDOWN: 662 case KSP_DIVERGED_INDEFINITE_MAT: 663 case KSP_DIVERGED_INDEFINITE_PC: 664 case KSP_CONVERGED_CG_NEG_CURVE: 665 /* Matrix or preconditioner is indefinite; increase perturbation */ 666 if (bnk->pert <= 0.0) { 667 /* Initialize the perturbation */ 668 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 669 if (bnk->is_gltr) { 670 ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 671 bnk->pert = PetscMax(bnk->pert, -e_min); 672 } 673 } else { 674 /* Increase the perturbation */ 675 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 676 } 677 break; 678 679 default: 680 /* Newton step computation is good; decrease perturbation */ 681 bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm); 682 if (bnk->pert < bnk->pmin) { 683 bnk->pert = 0.0; 684 } 685 break; 686 } 687 *stepType = BNK_NEWTON; 688 } 689 PetscFunctionReturn(0); 690 } 691 692 /*------------------------------------------------------------*/ 693 694 /* Routine for performing a bound-projected More-Thuente line search. 695 696 Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the 697 Newton step does not produce a valid step length. 698 */ 699 700 PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason) 701 { 702 TAO_BNK *bnk = (TAO_BNK *)tao->data; 703 PetscErrorCode ierr; 704 TaoLineSearchConvergedReason ls_reason; 705 706 PetscReal e_min, gdx, delta; 707 PetscInt bfgsUpdates; 708 709 PetscFunctionBegin; 710 /* Perform the linesearch */ 711 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 712 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 713 714 while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && (stepType != BNK_GRADIENT || stepType !=BNK_SCALED_GRADIENT)) { 715 /* Linesearch failed, revert solution */ 716 bnk->f = bnk->fold; 717 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 718 ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 719 720 switch(stepType) { 721 case BNK_NEWTON: 722 /* Failed to obtain acceptable iterate with Newton step 723 Update the perturbation for next time */ 724 if (bnk->pert <= 0.0) { 725 /* Initialize the perturbation */ 726 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 727 if (bnk->is_gltr) { 728 ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 729 bnk->pert = PetscMax(bnk->pert, -e_min); 730 } 731 } else { 732 /* Increase the perturbation */ 733 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 734 } 735 736 if (BNK_PC_BFGS != bnk->pc_type) { 737 /* We don't have the bfgs matrix around and being updated 738 Must use gradient direction in this case */ 739 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 740 stepType = BNK_GRADIENT; 741 } else { 742 /* Attempt to use the BFGS direction */ 743 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 744 /* Check for success (descent direction) 745 NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */ 746 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 747 if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 748 /* BFGS direction is not descent or direction produced not a number 749 We can assert bfgsUpdates > 1 in this case 750 Use steepest descent direction (scaled) */ 751 delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 752 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 753 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 754 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 755 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 756 757 bfgsUpdates = 1; 758 stepType = BNK_SCALED_GRADIENT; 759 } else { 760 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 761 if (1 == bfgsUpdates) { 762 /* The first BFGS direction is always the scaled gradient */ 763 stepType = BNK_SCALED_GRADIENT; 764 } else { 765 stepType = BNK_BFGS; 766 } 767 } 768 } 769 break; 770 771 case BNK_BFGS: 772 /* Can only enter if pc_type == BNK_PC_BFGS 773 Failed to obtain acceptable iterate with BFGS step 774 Attempt to use the scaled gradient direction */ 775 delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm); 776 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 777 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 778 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 779 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 780 781 bfgsUpdates = 1; 782 stepType = BNK_SCALED_GRADIENT; 783 break; 784 785 case BNK_SCALED_GRADIENT: 786 /* Can only enter if pc_type == BNK_PC_BFGS 787 The scaled gradient step did not produce a new iterate; 788 attemp to use the gradient direction. 789 Need to make sure we are not using a different diagonal scaling */ 790 ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 791 ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 792 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 793 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 794 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 795 796 bfgsUpdates = 1; 797 stepType = BNK_GRADIENT; 798 break; 799 } 800 /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 801 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 802 ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 803 804 /* Perform one last line search with the fall-back step */ 805 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr); 806 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 807 } 808 *reason = ls_reason; 809 PetscFunctionReturn(0); 810 } 811 812 /*------------------------------------------------------------*/ 813 814 /* Routine for updating the trust radius. 815 816 Function features three different update methods: 817 1) Line-search step length based 818 2) Predicted decrease on the CG quadratic model 819 3) Interpolation 820 */ 821 822 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept) 823 { 824 TAO_BNK *bnk = (TAO_BNK *)tao->data; 825 PetscErrorCode ierr; 826 827 PetscReal step, kappa; 828 PetscReal gdx, tau_1, tau_2, tau_min, tau_max; 829 830 PetscFunctionBegin; 831 /* Update trust region radius */ 832 *accept = PETSC_FALSE; 833 switch(updateType) { 834 case BNK_UPDATE_STEP: 835 *accept = PETSC_TRUE; /* always accept here because line search succeeded */ 836 if (stepType == BNK_NEWTON) { 837 ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr); 838 if (step < bnk->nu1) { 839 /* Very bad step taken; reduce radius */ 840 tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 841 } else if (step < bnk->nu2) { 842 /* Reasonably bad step taken; reduce radius */ 843 tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust); 844 } else if (step < bnk->nu3) { 845 /* Reasonable step was taken; leave radius alone */ 846 if (bnk->omega3 < 1.0) { 847 tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust); 848 } else if (bnk->omega3 > 1.0) { 849 tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust); 850 } 851 } else if (step < bnk->nu4) { 852 /* Full step taken; increase the radius */ 853 tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust); 854 } else { 855 /* More than full step taken; increase the radius */ 856 tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust); 857 } 858 } else { 859 /* Newton step was not good; reduce the radius */ 860 tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust); 861 } 862 break; 863 864 case BNK_UPDATE_REDUCTION: 865 if (stepType == BNK_NEWTON) { 866 if (prered < 0.0) { 867 /* The predicted reduction has the wrong sign. This cannot 868 happen in infinite precision arithmetic. Step should 869 be rejected! */ 870 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 871 } 872 else { 873 if (PetscIsInfOrNanReal(actred)) { 874 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 875 } else { 876 if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && 877 (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) { 878 kappa = 1.0; 879 } 880 else { 881 kappa = actred / prered; 882 } 883 884 /* Accept or reject the step and update radius */ 885 if (kappa < bnk->eta1) { 886 /* Reject the step */ 887 tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm); 888 } 889 else { 890 /* Accept the step */ 891 *accept = PETSC_TRUE; 892 /* Update the trust region radius only if the computed step is at the trust radius boundary */ 893 if (bnk->dnorm == tao->trust) { 894 if (kappa < bnk->eta2) { 895 /* Marginal bad step */ 896 tao->trust = bnk->alpha2 * tao->trust; 897 } 898 else if (kappa < bnk->eta3) { 899 /* Reasonable step */ 900 tao->trust = bnk->alpha3 * tao->trust; 901 } 902 else if (kappa < bnk->eta4) { 903 /* Good step */ 904 tao->trust = bnk->alpha4 * tao->trust; 905 } 906 else { 907 /* Very good step */ 908 tao->trust = bnk->alpha5 * tao->trust; 909 } 910 } 911 } 912 } 913 } 914 } else { 915 /* Newton step was not good; reduce the radius */ 916 tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust); 917 } 918 break; 919 920 default: 921 if (stepType == BNK_NEWTON) { 922 if (prered < 0.0) { 923 /* The predicted reduction has the wrong sign. This cannot */ 924 /* happen in infinite precision arithmetic. Step should */ 925 /* be rejected! */ 926 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 927 } else { 928 if (PetscIsInfOrNanReal(actred)) { 929 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 930 } else { 931 if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) { 932 kappa = 1.0; 933 } else { 934 kappa = actred / prered; 935 } 936 937 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 938 tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred); 939 tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred); 940 tau_min = PetscMin(tau_1, tau_2); 941 tau_max = PetscMax(tau_1, tau_2); 942 943 if (kappa >= 1.0 - bnk->mu1) { 944 /* Great agreement */ 945 *accept = PETSC_TRUE; 946 if (tau_max < 1.0) { 947 tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 948 } else if (tau_max > bnk->gamma4) { 949 tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm); 950 } else { 951 tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 952 } 953 } else if (kappa >= 1.0 - bnk->mu2) { 954 /* Good agreement */ 955 *accept = PETSC_TRUE; 956 if (tau_max < bnk->gamma2) { 957 tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 958 } else if (tau_max > bnk->gamma3) { 959 tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm); 960 } else if (tau_max < 1.0) { 961 tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 962 } else { 963 tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm); 964 } 965 } else { 966 /* Not good agreement */ 967 if (tau_min > 1.0) { 968 tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm); 969 } else if (tau_max < bnk->gamma1) { 970 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 971 } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) { 972 tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm); 973 } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) { 974 tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm); 975 } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) { 976 tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm); 977 } else { 978 tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm); 979 } 980 } 981 } 982 } 983 } else { 984 /* Newton step was not good; reduce the radius */ 985 tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust); 986 } 987 break; 988 } 989 /* Make sure the radius does not violate min and max settings */ 990 tao->trust = PetscMin(tao->trust, bnk->max_radius); 991 tao->trust = PetscMax(tao->trust, bnk->min_radius); 992 PetscFunctionReturn(0); 993 } 994 995 /* ---------------------------------------------------------- */ 996 997 PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType) 998 { 999 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1000 1001 PetscFunctionBegin; 1002 switch (stepType) { 1003 case BNK_NEWTON: 1004 ++bnk->newt; 1005 break; 1006 case BNK_BFGS: 1007 ++bnk->bfgs; 1008 break; 1009 case BNK_SCALED_GRADIENT: 1010 ++bnk->sgrad; 1011 break; 1012 case BNK_GRADIENT: 1013 ++bnk->grad; 1014 break; 1015 default: 1016 break; 1017 } 1018 PetscFunctionReturn(0); 1019 } 1020 1021 /* ---------------------------------------------------------- */ 1022 1023 PetscErrorCode TaoSetUp_BNK(Tao tao) 1024 { 1025 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1026 PetscErrorCode ierr; 1027 KSPType ksp_type; 1028 1029 PetscFunctionBegin; 1030 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 1031 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 1032 if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);} 1033 if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);} 1034 if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);} 1035 if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);} 1036 if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);} 1037 if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);} 1038 if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);} 1039 if (!bnk->G_inactive) {ierr = VecDuplicate(tao->solution,&bnk->G_inactive);CHKERRQ(ierr);} 1040 if (!bnk->Diag_min) {ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);} 1041 if (!bnk->Diag_max) {ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);} 1042 bnk->Diag = 0; 1043 bnk->inactive_work = 0; 1044 bnk->active_work = 0; 1045 bnk->inactive_idx = 0; 1046 bnk->active_idx = 0; 1047 bnk->active_lower = 0; 1048 bnk->active_upper = 0; 1049 bnk->active_fixed = 0; 1050 bnk->M = 0; 1051 bnk->H_inactive = 0; 1052 bnk->Hpre_inactive = 0; 1053 ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 1054 ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 1055 ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 1056 ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 1060 /*------------------------------------------------------------*/ 1061 1062 static PetscErrorCode TaoDestroy_BNK(Tao tao) 1063 { 1064 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1065 PetscErrorCode ierr; 1066 1067 PetscFunctionBegin; 1068 if (tao->setupcalled) { 1069 ierr = VecDestroy(&bnk->W);CHKERRQ(ierr); 1070 ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr); 1071 ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr); 1072 ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr); 1073 ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr); 1074 ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr); 1075 ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr); 1076 ierr = VecDestroy(&bnk->G_inactive);CHKERRQ(ierr); 1077 ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr); 1078 ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr); 1079 } 1080 ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr); 1081 ierr = MatDestroy(&bnk->M);CHKERRQ(ierr); 1082 if (bnk->Hpre_inactive != bnk->H_inactive) {ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);} 1083 ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 1084 ierr = PetscFree(tao->data);CHKERRQ(ierr); 1085 PetscFunctionReturn(0); 1086 } 1087 1088 /*------------------------------------------------------------*/ 1089 1090 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao) 1091 { 1092 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1093 PetscErrorCode ierr; 1094 1095 PetscFunctionBegin; 1096 ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr); 1097 ierr = PetscOptionsEList("-tao_bnk_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr); 1098 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); 1099 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); 1100 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); 1101 ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr); 1102 ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 1103 ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 1104 ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 1105 ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 1106 ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 1107 ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 1108 ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 1109 ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 1110 ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 1111 ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 1112 ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 1113 ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 1114 ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 1115 ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 1116 ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 1117 ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 1118 ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 1119 ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 1120 ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 1121 ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 1122 ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 1123 ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 1124 ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 1125 ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 1126 ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 1127 ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 1128 ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 1129 ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 1130 ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr); 1131 ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr); 1132 ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr); 1133 ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr); 1134 ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr); 1135 ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr); 1136 ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr); 1137 ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 1138 ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 1139 ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 1140 ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 1141 ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 1142 ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 1143 ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 1144 ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 1145 ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 1146 ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 1147 ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr); 1148 ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr); 1149 ierr = PetscOptionsTail();CHKERRQ(ierr); 1150 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1151 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1152 PetscFunctionReturn(0); 1153 } 1154 1155 /*------------------------------------------------------------*/ 1156 1157 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer) 1158 { 1159 TAO_BNK *bnk = (TAO_BNK *)tao->data; 1160 PetscInt nrejects; 1161 PetscBool isascii; 1162 PetscErrorCode ierr; 1163 1164 PetscFunctionBegin; 1165 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1166 if (isascii) { 1167 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1168 if (BNK_PC_BFGS == bnk->pc_type && bnk->M) { 1169 ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr); 1170 ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1171 } 1172 ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr); 1173 ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr); 1174 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr); 1175 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr); 1176 ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr); 1177 ierr = PetscViewerASCIIPrintf(viewer, " atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr); 1178 ierr = PetscViewerASCIIPrintf(viewer, " rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr); 1179 ierr = PetscViewerASCIIPrintf(viewer, " ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr); 1180 ierr = PetscViewerASCIIPrintf(viewer, " negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr); 1181 ierr = PetscViewerASCIIPrintf(viewer, " dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr); 1182 ierr = PetscViewerASCIIPrintf(viewer, " iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr); 1183 ierr = PetscViewerASCIIPrintf(viewer, " othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr); 1184 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1185 } 1186 PetscFunctionReturn(0); 1187 } 1188 1189 /* ---------------------------------------------------------- */ 1190 1191 /*MC 1192 TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms. 1193 At each iteration, the BNK methods solve the symmetric 1194 system of equations to obtain the step diretion dk: 1195 Hk dk = -gk 1196 for free variables only. The step can be globalized either through 1197 trust-region methods, or a line search, or a heuristic mixture of both. 1198 1199 Options Database Keys: 1200 + -tao_bnk_pc_type - "none","ahess","bfgs","petsc" 1201 . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs" 1202 . -tao_bnk_init_type - "constant","direction","interpolation" 1203 . -tao_bnk_update_type - "step","direction","interpolation" 1204 . -tao_bnk_as_type - "none","bertsekas" 1205 . -tao_bnk_sval - perturbation starting value 1206 . -tao_bnk_imin - minimum initial perturbation 1207 . -tao_bnk_imax - maximum initial perturbation 1208 . -tao_bnk_pmin - minimum perturbation 1209 . -tao_bnk_pmax - maximum perturbation 1210 . -tao_bnk_pgfac - growth factor 1211 . -tao_bnk_psfac - shrink factor 1212 . -tao_bnk_imfac - initial merit factor 1213 . -tao_bnk_pmgfac - merit growth factor 1214 . -tao_bnk_pmsfac - merit shrink factor 1215 . -tao_bnk_eta1 - poor steplength; reduce radius 1216 . -tao_bnk_eta2 - reasonable steplength; leave radius 1217 . -tao_bnk_eta3 - good steplength; increase readius 1218 . -tao_bnk_eta4 - excellent steplength; greatly increase radius 1219 . -tao_bnk_alpha1 - alpha1 reduction 1220 . -tao_bnk_alpha2 - alpha2 reduction 1221 . -tao_bnk_alpha3 - alpha3 reduction 1222 . -tao_bnk_alpha4 - alpha4 reduction 1223 . -tao_bnk_alpha - alpha5 reduction 1224 . -tao_bnk_mu1 - mu1 interpolation update 1225 . -tao_bnk_mu2 - mu2 interpolation update 1226 . -tao_bnk_gamma1 - gamma1 interpolation update 1227 . -tao_bnk_gamma2 - gamma2 interpolation update 1228 . -tao_bnk_gamma3 - gamma3 interpolation update 1229 . -tao_bnk_gamma4 - gamma4 interpolation update 1230 . -tao_bnk_theta - theta interpolation update 1231 . -tao_bnk_omega1 - omega1 step update 1232 . -tao_bnk_omega2 - omega2 step update 1233 . -tao_bnk_omega3 - omega3 step update 1234 . -tao_bnk_omega4 - omega4 step update 1235 . -tao_bnk_omega5 - omega5 step update 1236 . -tao_bnk_mu1_i - mu1 interpolation init factor 1237 . -tao_bnk_mu2_i - mu2 interpolation init factor 1238 . -tao_bnk_gamma1_i - gamma1 interpolation init factor 1239 . -tao_bnk_gamma2_i - gamma2 interpolation init factor 1240 . -tao_bnk_gamma3_i - gamma3 interpolation init factor 1241 . -tao_bnk_gamma4_i - gamma4 interpolation init factor 1242 . -tao_bnk_theta_i - theta interpolation init factor 1243 - -tao_bnk_bound_tol - initial tolerance used in estimating bounded active variables 1244 1245 Level: beginner 1246 M*/ 1247 1248 PetscErrorCode TaoCreate_BNK(Tao tao) 1249 { 1250 TAO_BNK *bnk; 1251 const char *morethuente_type = TAOLINESEARCHMT; 1252 PetscErrorCode ierr; 1253 1254 PetscFunctionBegin; 1255 ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr); 1256 1257 tao->ops->setup = TaoSetUp_BNK; 1258 tao->ops->view = TaoView_BNK; 1259 tao->ops->setfromoptions = TaoSetFromOptions_BNK; 1260 tao->ops->destroy = TaoDestroy_BNK; 1261 1262 /* Override default settings (unless already changed) */ 1263 if (!tao->max_it_changed) tao->max_it = 50; 1264 if (!tao->trust0_changed) tao->trust0 = 100.0; 1265 1266 tao->data = (void*)bnk; 1267 1268 /* Hessian shifting parameters */ 1269 bnk->sval = 0.0; 1270 bnk->imin = 1.0e-4; 1271 bnk->imax = 1.0e+2; 1272 bnk->imfac = 1.0e-1; 1273 1274 bnk->pmin = 1.0e-12; 1275 bnk->pmax = 1.0e+2; 1276 bnk->pgfac = 1.0e+1; 1277 bnk->psfac = 4.0e-1; 1278 bnk->pmgfac = 1.0e-1; 1279 bnk->pmsfac = 1.0e-1; 1280 1281 /* Default values for trust-region radius update based on steplength */ 1282 bnk->nu1 = 0.25; 1283 bnk->nu2 = 0.50; 1284 bnk->nu3 = 1.00; 1285 bnk->nu4 = 1.25; 1286 1287 bnk->omega1 = 0.25; 1288 bnk->omega2 = 0.50; 1289 bnk->omega3 = 1.00; 1290 bnk->omega4 = 2.00; 1291 bnk->omega5 = 4.00; 1292 1293 /* Default values for trust-region radius update based on reduction */ 1294 bnk->eta1 = 1.0e-4; 1295 bnk->eta2 = 0.25; 1296 bnk->eta3 = 0.50; 1297 bnk->eta4 = 0.90; 1298 1299 bnk->alpha1 = 0.25; 1300 bnk->alpha2 = 0.50; 1301 bnk->alpha3 = 1.00; 1302 bnk->alpha4 = 2.00; 1303 bnk->alpha5 = 4.00; 1304 1305 /* Default values for trust-region radius update based on interpolation */ 1306 bnk->mu1 = 0.10; 1307 bnk->mu2 = 0.50; 1308 1309 bnk->gamma1 = 0.25; 1310 bnk->gamma2 = 0.50; 1311 bnk->gamma3 = 2.00; 1312 bnk->gamma4 = 4.00; 1313 1314 bnk->theta = 0.05; 1315 1316 /* Default values for trust region initialization based on interpolation */ 1317 bnk->mu1_i = 0.35; 1318 bnk->mu2_i = 0.50; 1319 1320 bnk->gamma1_i = 0.0625; 1321 bnk->gamma2_i = 0.5; 1322 bnk->gamma3_i = 2.0; 1323 bnk->gamma4_i = 5.0; 1324 1325 bnk->theta_i = 0.25; 1326 1327 /* Remaining parameters */ 1328 bnk->min_radius = 1.0e-10; 1329 bnk->max_radius = 1.0e10; 1330 bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 1331 bnk->as_tol = 1.0e-3; 1332 bnk->as_step = 1.0e-3; 1333 bnk->dmin = 1.0e-6; 1334 bnk->dmax = 1.0e6; 1335 1336 bnk->pc_type = BNK_PC_BFGS; 1337 bnk->bfgs_scale_type = BFGS_SCALE_PHESS; 1338 bnk->init_type = BNK_INIT_INTERPOLATION; 1339 bnk->update_type = BNK_UPDATE_INTERPOLATION; 1340 bnk->as_type = BNK_AS_BERTSEKAS; 1341 1342 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1343 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 1344 ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1345 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1346 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 1347 1348 /* Set linear solver to default for symmetric matrices */ 1349 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1350 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 1351 ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 1352 ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 1353 PetscFunctionReturn(0); 1354 } 1355