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