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