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