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