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