1 #include <petsctaolinesearch.h> 2 #include <../src/tao/matrix/lmvmmat.h> 3 #include <../src/tao/unconstrained/impls/nls/nls.h> 4 5 #include <petscksp.h> 6 #include <petscpc.h> 7 #include <petsc-private/kspimpl.h> 8 #include <petsc-private/pcimpl.h> 9 10 #define NLS_KSP_CG 0 11 #define NLS_KSP_NASH 1 12 #define NLS_KSP_STCG 2 13 #define NLS_KSP_GLTR 3 14 #define NLS_KSP_PETSC 4 15 #define NLS_KSP_TYPES 5 16 17 #define NLS_PC_NONE 0 18 #define NLS_PC_AHESS 1 19 #define NLS_PC_BFGS 2 20 #define NLS_PC_PETSC 3 21 #define NLS_PC_TYPES 4 22 23 #define BFGS_SCALE_AHESS 0 24 #define BFGS_SCALE_PHESS 1 25 #define BFGS_SCALE_BFGS 2 26 #define BFGS_SCALE_TYPES 3 27 28 #define NLS_INIT_CONSTANT 0 29 #define NLS_INIT_DIRECTION 1 30 #define NLS_INIT_INTERPOLATION 2 31 #define NLS_INIT_TYPES 3 32 33 #define NLS_UPDATE_STEP 0 34 #define NLS_UPDATE_REDUCTION 1 35 #define NLS_UPDATE_INTERPOLATION 2 36 #define NLS_UPDATE_TYPES 3 37 38 static const char *NLS_KSP[64] = {"cg", "nash", "stcg", "gltr", "petsc"}; 39 40 static const char *NLS_PC[64] = {"none", "ahess", "bfgs", "petsc"}; 41 42 static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"}; 43 44 static const char *NLS_INIT[64] = {"constant", "direction", "interpolation"}; 45 46 static const char *NLS_UPDATE[64] = {"step", "reduction", "interpolation"}; 47 48 static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x); 49 /* Routine for BFGS preconditioner 50 51 52 Implements Newton's Method with a line search approach for solving 53 unconstrained minimization problems. A More'-Thuente line search 54 is used to guarantee that the bfgs preconditioner remains positive 55 definite. 56 57 The method can shift the Hessian matrix. The shifting procedure is 58 adapted from the PATH algorithm for solving complementarity 59 problems. 60 61 The linear system solve should be done with a conjugate gradient 62 method, although any method can be used. */ 63 64 #define NLS_NEWTON 0 65 #define NLS_BFGS 1 66 #define NLS_SCALED_GRADIENT 2 67 #define NLS_GRADIENT 3 68 69 #undef __FUNCT__ 70 #define __FUNCT__ "TaoSolve_NLS" 71 static PetscErrorCode TaoSolve_NLS(Tao tao) 72 { 73 PetscErrorCode ierr; 74 TAO_NLS *nlsP = (TAO_NLS *)tao->data; 75 PC pc; 76 KSPConvergedReason ksp_reason; 77 TaoLineSearchConvergedReason ls_reason; 78 TaoConvergedReason reason; 79 80 PetscReal fmin, ftrial, f_full, prered, actred, kappa, sigma; 81 PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius; 82 PetscReal f, fold, gdx, gnorm, pert; 83 PetscReal step = 1.0; 84 PetscReal delta; 85 PetscReal norm_d = 0.0, e_min; 86 87 PetscInt stepType; 88 PetscInt iter = 0; 89 PetscInt bfgsUpdates = 0; 90 PetscInt n,N,kspits; 91 PetscInt needH; 92 93 PetscInt i_max = 5; 94 PetscInt j_max = 1; 95 PetscInt i, j; 96 97 PetscFunctionBegin; 98 if (tao->XL || tao->XU || tao->ops->computebounds) { 99 ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr); 100 } 101 102 /* Initialized variables */ 103 pert = nlsP->sval; 104 105 /* Number of times ksp stopped because of these reasons */ 106 nlsP->ksp_atol = 0; 107 nlsP->ksp_rtol = 0; 108 nlsP->ksp_dtol = 0; 109 nlsP->ksp_ctol = 0; 110 nlsP->ksp_negc = 0; 111 nlsP->ksp_iter = 0; 112 nlsP->ksp_othr = 0; 113 114 /* Modify the linear solver to a trust region method if desired */ 115 switch(nlsP->ksp_type) { 116 case NLS_KSP_CG: 117 ierr = KSPSetType(tao->ksp, KSPCG);CHKERRQ(ierr); 118 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 119 break; 120 121 case NLS_KSP_NASH: 122 ierr = KSPSetType(tao->ksp, KSPNASH);CHKERRQ(ierr); 123 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 124 break; 125 126 case NLS_KSP_STCG: 127 ierr = KSPSetType(tao->ksp, KSPSTCG);CHKERRQ(ierr); 128 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 129 break; 130 131 case NLS_KSP_GLTR: 132 ierr = KSPSetType(tao->ksp, KSPGLTR);CHKERRQ(ierr); 133 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 134 break; 135 136 default: 137 /* Use the method set by the ksp_type */ 138 break; 139 } 140 141 /* Initialize trust-region radius when using nash, stcg, or gltr 142 Will be reset during the first iteration */ 143 if (NLS_KSP_NASH == nlsP->ksp_type) { 144 ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 145 } else if (NLS_KSP_STCG == nlsP->ksp_type) { 146 ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 147 } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 148 ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 149 } 150 151 if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 152 tao->trust = tao->trust0; 153 154 if (tao->trust < 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial radius negative"); 155 156 /* Modify the radius if it is too large or small */ 157 tao->trust = PetscMax(tao->trust, nlsP->min_radius); 158 tao->trust = PetscMin(tao->trust, nlsP->max_radius); 159 } 160 161 /* Get vectors we will need */ 162 if (NLS_PC_BFGS == nlsP->pc_type && !nlsP->M) { 163 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 164 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 165 ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&nlsP->M);CHKERRQ(ierr); 166 ierr = MatLMVMAllocateVectors(nlsP->M,tao->solution);CHKERRQ(ierr); 167 } 168 169 /* Check convergence criteria */ 170 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 171 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 172 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 173 needH = 1; 174 175 ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 176 if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 177 178 /* create vectors for the limited memory preconditioner */ 179 if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_BFGS != nlsP->bfgs_scale_type)) { 180 if (!nlsP->Diag) { 181 ierr = VecDuplicate(tao->solution,&nlsP->Diag);CHKERRQ(ierr); 182 } 183 } 184 185 /* Modify the preconditioner to use the bfgs approximation */ 186 ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr); 187 switch(nlsP->pc_type) { 188 case NLS_PC_NONE: 189 ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr); 190 if (pc->ops->setfromoptions) { 191 (*pc->ops->setfromoptions)(pc); 192 } 193 break; 194 195 case NLS_PC_AHESS: 196 ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr); 197 if (pc->ops->setfromoptions) { 198 (*pc->ops->setfromoptions)(pc); 199 } 200 ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr); 201 break; 202 203 case NLS_PC_BFGS: 204 ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr); 205 if (pc->ops->setfromoptions) { 206 (*pc->ops->setfromoptions)(pc); 207 } 208 ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr); 209 ierr = PCShellSetContext(pc, nlsP->M);CHKERRQ(ierr); 210 ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr); 211 break; 212 213 default: 214 /* Use the pc method set by pc_type */ 215 break; 216 } 217 218 /* Initialize trust-region radius. The initialization is only performed 219 when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */ 220 if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 221 switch(nlsP->init_type) { 222 case NLS_INIT_CONSTANT: 223 /* Use the initial radius specified */ 224 break; 225 226 case NLS_INIT_INTERPOLATION: 227 /* Use the initial radius specified */ 228 max_radius = 0.0; 229 230 for (j = 0; j < j_max; ++j) { 231 fmin = f; 232 sigma = 0.0; 233 234 if (needH) { 235 ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 236 needH = 0; 237 } 238 239 for (i = 0; i < i_max; ++i) { 240 ierr = VecCopy(tao->solution,nlsP->W);CHKERRQ(ierr); 241 ierr = VecAXPY(nlsP->W,-tao->trust/gnorm,tao->gradient);CHKERRQ(ierr); 242 ierr = TaoComputeObjective(tao, nlsP->W, &ftrial);CHKERRQ(ierr); 243 if (PetscIsInfOrNanReal(ftrial)) { 244 tau = nlsP->gamma1_i; 245 } else { 246 if (ftrial < fmin) { 247 fmin = ftrial; 248 sigma = -tao->trust / gnorm; 249 } 250 251 ierr = MatMult(tao->hessian, tao->gradient, nlsP->D);CHKERRQ(ierr); 252 ierr = VecDot(tao->gradient, nlsP->D, &prered);CHKERRQ(ierr); 253 254 prered = tao->trust * (gnorm - 0.5 * tao->trust * prered / (gnorm * gnorm)); 255 actred = f - ftrial; 256 if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 257 kappa = 1.0; 258 } else { 259 kappa = actred / prered; 260 } 261 262 tau_1 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust + (1.0 - nlsP->theta_i) * prered - actred); 263 tau_2 = nlsP->theta_i * gnorm * tao->trust / (nlsP->theta_i * gnorm * tao->trust - (1.0 + nlsP->theta_i) * prered + actred); 264 tau_min = PetscMin(tau_1, tau_2); 265 tau_max = PetscMax(tau_1, tau_2); 266 267 if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu1_i) { 268 /* Great agreement */ 269 max_radius = PetscMax(max_radius, tao->trust); 270 271 if (tau_max < 1.0) { 272 tau = nlsP->gamma3_i; 273 } else if (tau_max > nlsP->gamma4_i) { 274 tau = nlsP->gamma4_i; 275 } else if (tau_1 >= 1.0 && tau_1 <= nlsP->gamma4_i && tau_2 < 1.0) { 276 tau = tau_1; 277 } else if (tau_2 >= 1.0 && tau_2 <= nlsP->gamma4_i && tau_1 < 1.0) { 278 tau = tau_2; 279 } else { 280 tau = tau_max; 281 } 282 } else if (PetscAbsScalar(kappa - 1.0) <= nlsP->mu2_i) { 283 /* Good agreement */ 284 max_radius = PetscMax(max_radius, tao->trust); 285 286 if (tau_max < nlsP->gamma2_i) { 287 tau = nlsP->gamma2_i; 288 } else if (tau_max > nlsP->gamma3_i) { 289 tau = nlsP->gamma3_i; 290 } else { 291 tau = tau_max; 292 } 293 } else { 294 /* Not good agreement */ 295 if (tau_min > 1.0) { 296 tau = nlsP->gamma2_i; 297 } else if (tau_max < nlsP->gamma1_i) { 298 tau = nlsP->gamma1_i; 299 } else if ((tau_min < nlsP->gamma1_i) && (tau_max >= 1.0)) { 300 tau = nlsP->gamma1_i; 301 } else if ((tau_1 >= nlsP->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 302 tau = tau_1; 303 } else if ((tau_2 >= nlsP->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1_i) || (tau_2 >= 1.0))) { 304 tau = tau_2; 305 } else { 306 tau = tau_max; 307 } 308 } 309 } 310 tao->trust = tau * tao->trust; 311 } 312 313 if (fmin < f) { 314 f = fmin; 315 ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr); 316 ierr = TaoComputeGradient(tao,tao->solution,tao->gradient);CHKERRQ(ierr); 317 318 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 319 if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN"); 320 needH = 1; 321 322 ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, 1.0, &reason);CHKERRQ(ierr); 323 if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 324 } 325 } 326 tao->trust = PetscMax(tao->trust, max_radius); 327 328 /* Modify the radius if it is too large or small */ 329 tao->trust = PetscMax(tao->trust, nlsP->min_radius); 330 tao->trust = PetscMin(tao->trust, nlsP->max_radius); 331 break; 332 333 default: 334 /* Norm of the first direction will initialize radius */ 335 tao->trust = 0.0; 336 break; 337 } 338 } 339 340 /* Set initial scaling for the BFGS preconditioner 341 This step is done after computing the initial trust-region radius 342 since the function value may have decreased */ 343 if (NLS_PC_BFGS == nlsP->pc_type) { 344 if (f != 0.0) { 345 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 346 } else { 347 delta = 2.0 / (gnorm*gnorm); 348 } 349 ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr); 350 } 351 352 /* Set counter for gradient/reset steps*/ 353 nlsP->newt = 0; 354 nlsP->bfgs = 0; 355 nlsP->sgrad = 0; 356 nlsP->grad = 0; 357 358 /* Have not converged; continue with Newton method */ 359 while (reason == TAO_CONTINUE_ITERATING) { 360 ++iter; 361 tao->ksp_its=0; 362 363 /* Compute the Hessian */ 364 if (needH) { 365 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 366 needH = 0; 367 } 368 369 if ((NLS_PC_BFGS == nlsP->pc_type) && (BFGS_SCALE_AHESS == nlsP->bfgs_scale_type)) { 370 /* Obtain diagonal for the bfgs preconditioner */ 371 ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr); 372 ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr); 373 ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr); 374 ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr); 375 } 376 377 /* Shift the Hessian matrix */ 378 if (pert > 0) { 379 ierr = MatShift(tao->hessian, pert); 380 if (tao->hessian != tao->hessian_pre) { 381 ierr = MatShift(tao->hessian_pre, pert);CHKERRQ(ierr); 382 } 383 } 384 385 if (NLS_PC_BFGS == nlsP->pc_type) { 386 if (BFGS_SCALE_PHESS == nlsP->bfgs_scale_type) { 387 /* Obtain diagonal for the bfgs preconditioner */ 388 ierr = MatGetDiagonal(tao->hessian, nlsP->Diag);CHKERRQ(ierr); 389 ierr = VecAbs(nlsP->Diag);CHKERRQ(ierr); 390 ierr = VecReciprocal(nlsP->Diag);CHKERRQ(ierr); 391 ierr = MatLMVMSetScale(nlsP->M,nlsP->Diag);CHKERRQ(ierr); 392 } 393 /* Update the limited memory preconditioner */ 394 ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 395 ++bfgsUpdates; 396 } 397 398 /* Solve the Newton system of equations */ 399 ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 400 if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 401 402 if (NLS_KSP_NASH == nlsP->ksp_type) { 403 ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 404 } else if (NLS_KSP_STCG == nlsP->ksp_type) { 405 ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 406 } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 407 ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 408 } 409 410 ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 411 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 412 tao->ksp_its+=kspits; 413 tao->ksp_tot_its+=kspits; 414 415 if (NLS_KSP_NASH == nlsP->ksp_type) { 416 ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 417 } else if (NLS_KSP_STCG == nlsP->ksp_type) { 418 ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 419 } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 420 ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 421 } 422 423 if (0.0 == tao->trust) { 424 /* Radius was uninitialized; use the norm of the direction */ 425 if (norm_d > 0.0) { 426 tao->trust = norm_d; 427 428 /* Modify the radius if it is too large or small */ 429 tao->trust = PetscMax(tao->trust, nlsP->min_radius); 430 tao->trust = PetscMin(tao->trust, nlsP->max_radius); 431 } else { 432 /* The direction was bad; set radius to default value and re-solve 433 the trust-region subproblem to get a direction */ 434 tao->trust = tao->trust0; 435 436 /* Modify the radius if it is too large or small */ 437 tao->trust = PetscMax(tao->trust, nlsP->min_radius); 438 tao->trust = PetscMin(tao->trust, nlsP->max_radius); 439 440 if (NLS_KSP_NASH == nlsP->ksp_type) { 441 ierr = KSPNASHSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 442 } else if (NLS_KSP_STCG == nlsP->ksp_type) { 443 ierr = KSPSTCGSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 444 } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 445 ierr = KSPGLTRSetRadius(tao->ksp,nlsP->max_radius);CHKERRQ(ierr); 446 } 447 448 ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 449 ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr); 450 tao->ksp_its+=kspits; 451 tao->ksp_tot_its+=kspits; 452 if (NLS_KSP_NASH == nlsP->ksp_type) { 453 ierr = KSPNASHGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 454 } else if (NLS_KSP_STCG == nlsP->ksp_type) { 455 ierr = KSPSTCGGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 456 } else if (NLS_KSP_GLTR == nlsP->ksp_type) { 457 ierr = KSPGLTRGetNormD(tao->ksp,&norm_d);CHKERRQ(ierr); 458 } 459 460 if (norm_d == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero"); 461 } 462 } 463 } else { 464 ierr = KSPSolve(tao->ksp, tao->gradient, nlsP->D);CHKERRQ(ierr); 465 ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr); 466 tao->ksp_its += kspits; 467 tao->ksp_tot_its+=kspits; 468 } 469 ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 470 ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr); 471 if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (NLS_PC_BFGS == nlsP->pc_type) && (bfgsUpdates > 1)) { 472 /* Preconditioner is numerically indefinite; reset the 473 approximate if using BFGS preconditioning. */ 474 475 if (f != 0.0) { 476 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 477 } else { 478 delta = 2.0 / (gnorm*gnorm); 479 } 480 ierr = MatLMVMSetDelta(nlsP->M,delta);CHKERRQ(ierr); 481 ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 482 ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 483 bfgsUpdates = 1; 484 } 485 486 if (KSP_CONVERGED_ATOL == ksp_reason) { 487 ++nlsP->ksp_atol; 488 } else if (KSP_CONVERGED_RTOL == ksp_reason) { 489 ++nlsP->ksp_rtol; 490 } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) { 491 ++nlsP->ksp_ctol; 492 } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) { 493 ++nlsP->ksp_negc; 494 } else if (KSP_DIVERGED_DTOL == ksp_reason) { 495 ++nlsP->ksp_dtol; 496 } else if (KSP_DIVERGED_ITS == ksp_reason) { 497 ++nlsP->ksp_iter; 498 } else { 499 ++nlsP->ksp_othr; 500 } 501 502 /* Check for success (descent direction) */ 503 ierr = VecDot(nlsP->D, tao->gradient, &gdx);CHKERRQ(ierr); 504 if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 505 /* Newton step is not descent or direction produced Inf or NaN 506 Update the perturbation for next time */ 507 if (pert <= 0.0) { 508 /* Initialize the perturbation */ 509 pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 510 if (NLS_KSP_GLTR == nlsP->ksp_type) { 511 ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 512 pert = PetscMax(pert, -e_min); 513 } 514 } else { 515 /* Increase the perturbation */ 516 pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 517 } 518 519 if (NLS_PC_BFGS != nlsP->pc_type) { 520 /* We don't have the bfgs matrix around and updated 521 Must use gradient direction in this case */ 522 ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 523 ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 524 ++nlsP->grad; 525 stepType = NLS_GRADIENT; 526 } else { 527 /* Attempt to use the BFGS direction */ 528 ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 529 ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 530 531 /* Check for success (descent direction) */ 532 ierr = VecDot(tao->gradient, nlsP->D, &gdx);CHKERRQ(ierr); 533 if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 534 /* BFGS direction is not descent or direction produced not a number 535 We can assert bfgsUpdates > 1 in this case because 536 the first solve produces the scaled gradient direction, 537 which is guaranteed to be descent */ 538 539 /* Use steepest descent direction (scaled) */ 540 541 if (f != 0.0) { 542 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 543 } else { 544 delta = 2.0 / (gnorm*gnorm); 545 } 546 ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 547 ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 548 ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 549 ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 550 ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 551 552 bfgsUpdates = 1; 553 ++nlsP->sgrad; 554 stepType = NLS_SCALED_GRADIENT; 555 } else { 556 if (1 == bfgsUpdates) { 557 /* The first BFGS direction is always the scaled gradient */ 558 ++nlsP->sgrad; 559 stepType = NLS_SCALED_GRADIENT; 560 } else { 561 ++nlsP->bfgs; 562 stepType = NLS_BFGS; 563 } 564 } 565 } 566 } else { 567 /* Computed Newton step is descent */ 568 switch (ksp_reason) { 569 case KSP_DIVERGED_NANORINF: 570 case KSP_DIVERGED_BREAKDOWN: 571 case KSP_DIVERGED_INDEFINITE_MAT: 572 case KSP_DIVERGED_INDEFINITE_PC: 573 case KSP_CONVERGED_CG_NEG_CURVE: 574 /* Matrix or preconditioner is indefinite; increase perturbation */ 575 if (pert <= 0.0) { 576 /* Initialize the perturbation */ 577 pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 578 if (NLS_KSP_GLTR == nlsP->ksp_type) { 579 ierr = KSPGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr); 580 pert = PetscMax(pert, -e_min); 581 } 582 } else { 583 /* Increase the perturbation */ 584 pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 585 } 586 break; 587 588 default: 589 /* Newton step computation is good; decrease perturbation */ 590 pert = PetscMin(nlsP->psfac * pert, nlsP->pmsfac * gnorm); 591 if (pert < nlsP->pmin) { 592 pert = 0.0; 593 } 594 break; 595 } 596 597 ++nlsP->newt; 598 stepType = NLS_NEWTON; 599 } 600 601 /* Perform the linesearch */ 602 fold = f; 603 ierr = VecCopy(tao->solution, nlsP->Xold);CHKERRQ(ierr); 604 ierr = VecCopy(tao->gradient, nlsP->Gold);CHKERRQ(ierr); 605 606 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 607 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 608 609 while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != NLS_GRADIENT) { /* Linesearch failed */ 610 f = fold; 611 ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 612 ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 613 614 switch(stepType) { 615 case NLS_NEWTON: 616 /* Failed to obtain acceptable iterate with Newton 1step 617 Update the perturbation for next time */ 618 if (pert <= 0.0) { 619 /* Initialize the perturbation */ 620 pert = PetscMin(nlsP->imax, PetscMax(nlsP->imin, nlsP->imfac * gnorm)); 621 if (NLS_KSP_GLTR == nlsP->ksp_type) { 622 ierr = KSPGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 623 pert = PetscMax(pert, -e_min); 624 } 625 } else { 626 /* Increase the perturbation */ 627 pert = PetscMin(nlsP->pmax, PetscMax(nlsP->pgfac * pert, nlsP->pmgfac * gnorm)); 628 } 629 630 if (NLS_PC_BFGS != nlsP->pc_type) { 631 /* We don't have the bfgs matrix around and being updated 632 Must use gradient direction in this case */ 633 ierr = VecCopy(tao->gradient, nlsP->D);CHKERRQ(ierr); 634 ++nlsP->grad; 635 stepType = NLS_GRADIENT; 636 } else { 637 /* Attempt to use the BFGS direction */ 638 ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 639 /* Check for success (descent direction) */ 640 ierr = VecDot(tao->solution, nlsP->D, &gdx);CHKERRQ(ierr); 641 if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 642 /* BFGS direction is not descent or direction produced not a number 643 We can assert bfgsUpdates > 1 in this case 644 Use steepest descent direction (scaled) */ 645 646 if (f != 0.0) { 647 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 648 } else { 649 delta = 2.0 / (gnorm*gnorm); 650 } 651 ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 652 ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 653 ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 654 ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 655 656 bfgsUpdates = 1; 657 ++nlsP->sgrad; 658 stepType = NLS_SCALED_GRADIENT; 659 } else { 660 if (1 == bfgsUpdates) { 661 /* The first BFGS direction is always the scaled gradient */ 662 ++nlsP->sgrad; 663 stepType = NLS_SCALED_GRADIENT; 664 } else { 665 ++nlsP->bfgs; 666 stepType = NLS_BFGS; 667 } 668 } 669 } 670 break; 671 672 case NLS_BFGS: 673 /* Can only enter if pc_type == NLS_PC_BFGS 674 Failed to obtain acceptable iterate with BFGS step 675 Attempt to use the scaled gradient direction */ 676 677 if (f != 0.0) { 678 delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm); 679 } else { 680 delta = 2.0 / (gnorm*gnorm); 681 } 682 ierr = MatLMVMSetDelta(nlsP->M, delta);CHKERRQ(ierr); 683 ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 684 ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 685 ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 686 687 bfgsUpdates = 1; 688 ++nlsP->sgrad; 689 stepType = NLS_SCALED_GRADIENT; 690 break; 691 692 case NLS_SCALED_GRADIENT: 693 /* Can only enter if pc_type == NLS_PC_BFGS 694 The scaled gradient step did not produce a new iterate; 695 attemp to use the gradient direction. 696 Need to make sure we are not using a different diagonal scaling */ 697 698 ierr = MatLMVMSetScale(nlsP->M,0);CHKERRQ(ierr); 699 ierr = MatLMVMSetDelta(nlsP->M,1.0);CHKERRQ(ierr); 700 ierr = MatLMVMReset(nlsP->M);CHKERRQ(ierr); 701 ierr = MatLMVMUpdate(nlsP->M, tao->solution, tao->gradient);CHKERRQ(ierr); 702 ierr = MatLMVMSolve(nlsP->M, tao->gradient, nlsP->D);CHKERRQ(ierr); 703 704 bfgsUpdates = 1; 705 ++nlsP->grad; 706 stepType = NLS_GRADIENT; 707 break; 708 } 709 ierr = VecScale(nlsP->D, -1.0);CHKERRQ(ierr); 710 711 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, nlsP->D, &step, &ls_reason);CHKERRQ(ierr); 712 ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 713 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 714 } 715 716 if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 717 /* Failed to find an improving point */ 718 f = fold; 719 ierr = VecCopy(nlsP->Xold, tao->solution);CHKERRQ(ierr); 720 ierr = VecCopy(nlsP->Gold, tao->gradient);CHKERRQ(ierr); 721 step = 0.0; 722 reason = TAO_DIVERGED_LS_FAILURE; 723 tao->reason = TAO_DIVERGED_LS_FAILURE; 724 break; 725 } 726 727 /* Update trust region radius */ 728 if (NLS_KSP_NASH == nlsP->ksp_type || NLS_KSP_STCG == nlsP->ksp_type || NLS_KSP_GLTR == nlsP->ksp_type) { 729 switch(nlsP->update_type) { 730 case NLS_UPDATE_STEP: 731 if (stepType == NLS_NEWTON) { 732 if (step < nlsP->nu1) { 733 /* Very bad step taken; reduce radius */ 734 tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 735 } else if (step < nlsP->nu2) { 736 /* Reasonably bad step taken; reduce radius */ 737 tao->trust = nlsP->omega2 * PetscMin(norm_d, tao->trust); 738 } else if (step < nlsP->nu3) { 739 /* Reasonable step was taken; leave radius alone */ 740 if (nlsP->omega3 < 1.0) { 741 tao->trust = nlsP->omega3 * PetscMin(norm_d, tao->trust); 742 } else if (nlsP->omega3 > 1.0) { 743 tao->trust = PetscMax(nlsP->omega3 * norm_d, tao->trust); 744 } 745 } else if (step < nlsP->nu4) { 746 /* Full step taken; increase the radius */ 747 tao->trust = PetscMax(nlsP->omega4 * norm_d, tao->trust); 748 } else { 749 /* More than full step taken; increase the radius */ 750 tao->trust = PetscMax(nlsP->omega5 * norm_d, tao->trust); 751 } 752 } else { 753 /* Newton step was not good; reduce the radius */ 754 tao->trust = nlsP->omega1 * PetscMin(norm_d, tao->trust); 755 } 756 break; 757 758 case NLS_UPDATE_REDUCTION: 759 if (stepType == NLS_NEWTON) { 760 /* Get predicted reduction */ 761 762 if (NLS_KSP_STCG == nlsP->ksp_type) { 763 ierr = KSPSTCGGetObjFcn(tao->ksp,&prered); 764 } else if (NLS_KSP_NASH == nlsP->ksp_type) { 765 ierr = KSPNASHGetObjFcn(tao->ksp,&prered); 766 } else { 767 ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 768 } 769 770 if (prered >= 0.0) { 771 /* The predicted reduction has the wrong sign. This cannot */ 772 /* happen in infinite precision arithmetic. Step should */ 773 /* be rejected! */ 774 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 775 } else { 776 if (PetscIsInfOrNanReal(f_full)) { 777 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 778 } else { 779 /* Compute and actual reduction */ 780 actred = fold - f_full; 781 prered = -prered; 782 if ((PetscAbsScalar(actred) <= nlsP->epsilon) && 783 (PetscAbsScalar(prered) <= nlsP->epsilon)) { 784 kappa = 1.0; 785 } else { 786 kappa = actred / prered; 787 } 788 789 /* Accept of reject the step and update radius */ 790 if (kappa < nlsP->eta1) { 791 /* Very bad step */ 792 tao->trust = nlsP->alpha1 * PetscMin(tao->trust, norm_d); 793 } else if (kappa < nlsP->eta2) { 794 /* Marginal bad step */ 795 tao->trust = nlsP->alpha2 * PetscMin(tao->trust, norm_d); 796 } else if (kappa < nlsP->eta3) { 797 /* Reasonable step */ 798 if (nlsP->alpha3 < 1.0) { 799 tao->trust = nlsP->alpha3 * PetscMin(norm_d, tao->trust); 800 } else if (nlsP->alpha3 > 1.0) { 801 tao->trust = PetscMax(nlsP->alpha3 * norm_d, tao->trust); 802 } 803 } else if (kappa < nlsP->eta4) { 804 /* Good step */ 805 tao->trust = PetscMax(nlsP->alpha4 * norm_d, tao->trust); 806 } else { 807 /* Very good step */ 808 tao->trust = PetscMax(nlsP->alpha5 * norm_d, tao->trust); 809 } 810 } 811 } 812 } else { 813 /* Newton step was not good; reduce the radius */ 814 tao->trust = nlsP->alpha1 * PetscMin(norm_d, tao->trust); 815 } 816 break; 817 818 default: 819 if (stepType == NLS_NEWTON) { 820 821 if (NLS_KSP_STCG == nlsP->ksp_type) { 822 ierr = KSPSTCGGetObjFcn(tao->ksp,&prered); 823 } else if (NLS_KSP_NASH == nlsP->ksp_type) { 824 ierr = KSPNASHGetObjFcn(tao->ksp,&prered); 825 } else { 826 ierr = KSPGLTRGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr); 827 } 828 if (prered >= 0.0) { 829 /* The predicted reduction has the wrong sign. This cannot */ 830 /* happen in infinite precision arithmetic. Step should */ 831 /* be rejected! */ 832 tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 833 } else { 834 if (PetscIsInfOrNanReal(f_full)) { 835 tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 836 } else { 837 actred = fold - f_full; 838 prered = -prered; 839 if ((PetscAbsScalar(actred) <= nlsP->epsilon) && (PetscAbsScalar(prered) <= nlsP->epsilon)) { 840 kappa = 1.0; 841 } else { 842 kappa = actred / prered; 843 } 844 845 tau_1 = nlsP->theta * gdx / (nlsP->theta * gdx - (1.0 - nlsP->theta) * prered + actred); 846 tau_2 = nlsP->theta * gdx / (nlsP->theta * gdx + (1.0 + nlsP->theta) * prered - actred); 847 tau_min = PetscMin(tau_1, tau_2); 848 tau_max = PetscMax(tau_1, tau_2); 849 850 if (kappa >= 1.0 - nlsP->mu1) { 851 /* Great agreement */ 852 if (tau_max < 1.0) { 853 tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 854 } else if (tau_max > nlsP->gamma4) { 855 tao->trust = PetscMax(tao->trust, nlsP->gamma4 * norm_d); 856 } else { 857 tao->trust = PetscMax(tao->trust, tau_max * norm_d); 858 } 859 } else if (kappa >= 1.0 - nlsP->mu2) { 860 /* Good agreement */ 861 862 if (tau_max < nlsP->gamma2) { 863 tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 864 } else if (tau_max > nlsP->gamma3) { 865 tao->trust = PetscMax(tao->trust, nlsP->gamma3 * norm_d); 866 } else if (tau_max < 1.0) { 867 tao->trust = tau_max * PetscMin(tao->trust, norm_d); 868 } else { 869 tao->trust = PetscMax(tao->trust, tau_max * norm_d); 870 } 871 } else { 872 /* Not good agreement */ 873 if (tau_min > 1.0) { 874 tao->trust = nlsP->gamma2 * PetscMin(tao->trust, norm_d); 875 } else if (tau_max < nlsP->gamma1) { 876 tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 877 } else if ((tau_min < nlsP->gamma1) && (tau_max >= 1.0)) { 878 tao->trust = nlsP->gamma1 * PetscMin(tao->trust, norm_d); 879 } else if ((tau_1 >= nlsP->gamma1) && (tau_1 < 1.0) && ((tau_2 < nlsP->gamma1) || (tau_2 >= 1.0))) { 880 tao->trust = tau_1 * PetscMin(tao->trust, norm_d); 881 } else if ((tau_2 >= nlsP->gamma1) && (tau_2 < 1.0) && ((tau_1 < nlsP->gamma1) || (tau_2 >= 1.0))) { 882 tao->trust = tau_2 * PetscMin(tao->trust, norm_d); 883 } else { 884 tao->trust = tau_max * PetscMin(tao->trust, norm_d); 885 } 886 } 887 } 888 } 889 } else { 890 /* Newton step was not good; reduce the radius */ 891 tao->trust = nlsP->gamma1 * PetscMin(norm_d, tao->trust); 892 } 893 break; 894 } 895 896 /* The radius may have been increased; modify if it is too large */ 897 tao->trust = PetscMin(tao->trust, nlsP->max_radius); 898 } 899 900 /* Check for termination */ 901 ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr); 902 if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 903 needH = 1; 904 ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 905 } 906 PetscFunctionReturn(0); 907 } 908 909 /* ---------------------------------------------------------- */ 910 #undef __FUNCT__ 911 #define __FUNCT__ "TaoSetUp_NLS" 912 static PetscErrorCode TaoSetUp_NLS(Tao tao) 913 { 914 TAO_NLS *nlsP = (TAO_NLS *)tao->data; 915 PetscErrorCode ierr; 916 917 PetscFunctionBegin; 918 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 919 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 920 if (!nlsP->W) {ierr = VecDuplicate(tao->solution,&nlsP->W);CHKERRQ(ierr);} 921 if (!nlsP->D) {ierr = VecDuplicate(tao->solution,&nlsP->D);CHKERRQ(ierr);} 922 if (!nlsP->Xold) {ierr = VecDuplicate(tao->solution,&nlsP->Xold);CHKERRQ(ierr);} 923 if (!nlsP->Gold) {ierr = VecDuplicate(tao->solution,&nlsP->Gold);CHKERRQ(ierr);} 924 nlsP->Diag = 0; 925 nlsP->M = 0; 926 PetscFunctionReturn(0); 927 } 928 929 /*------------------------------------------------------------*/ 930 #undef __FUNCT__ 931 #define __FUNCT__ "TaoDestroy_NLS" 932 static PetscErrorCode TaoDestroy_NLS(Tao tao) 933 { 934 TAO_NLS *nlsP = (TAO_NLS *)tao->data; 935 PetscErrorCode ierr; 936 937 PetscFunctionBegin; 938 if (tao->setupcalled) { 939 ierr = VecDestroy(&nlsP->D);CHKERRQ(ierr); 940 ierr = VecDestroy(&nlsP->W);CHKERRQ(ierr); 941 ierr = VecDestroy(&nlsP->Xold);CHKERRQ(ierr); 942 ierr = VecDestroy(&nlsP->Gold);CHKERRQ(ierr); 943 } 944 ierr = VecDestroy(&nlsP->Diag);CHKERRQ(ierr); 945 ierr = MatDestroy(&nlsP->M);CHKERRQ(ierr); 946 ierr = PetscFree(tao->data);CHKERRQ(ierr); 947 PetscFunctionReturn(0); 948 } 949 950 /*------------------------------------------------------------*/ 951 #undef __FUNCT__ 952 #define __FUNCT__ "TaoSetFromOptions_NLS" 953 static PetscErrorCode TaoSetFromOptions_NLS(Tao tao) 954 { 955 TAO_NLS *nlsP = (TAO_NLS *)tao->data; 956 PetscErrorCode ierr; 957 958 PetscFunctionBegin; 959 ierr = PetscOptionsHead("Newton line search method for unconstrained optimization");CHKERRQ(ierr); 960 ierr = PetscOptionsEList("-tao_nls_ksp_type", "ksp type", "", NLS_KSP, NLS_KSP_TYPES, NLS_KSP[nlsP->ksp_type], &nlsP->ksp_type, 0);CHKERRQ(ierr); 961 ierr = PetscOptionsEList("-tao_nls_pc_type", "pc type", "", NLS_PC, NLS_PC_TYPES, NLS_PC[nlsP->pc_type], &nlsP->pc_type, 0);CHKERRQ(ierr); 962 ierr = PetscOptionsEList("-tao_nls_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[nlsP->bfgs_scale_type], &nlsP->bfgs_scale_type, 0);CHKERRQ(ierr); 963 ierr = PetscOptionsEList("-tao_nls_init_type", "radius initialization type", "", NLS_INIT, NLS_INIT_TYPES, NLS_INIT[nlsP->init_type], &nlsP->init_type, 0);CHKERRQ(ierr); 964 ierr = PetscOptionsEList("-tao_nls_update_type", "radius update type", "", NLS_UPDATE, NLS_UPDATE_TYPES, NLS_UPDATE[nlsP->update_type], &nlsP->update_type, 0);CHKERRQ(ierr); 965 ierr = PetscOptionsReal("-tao_nls_sval", "perturbation starting value", "", nlsP->sval, &nlsP->sval,NULL);CHKERRQ(ierr); 966 ierr = PetscOptionsReal("-tao_nls_imin", "minimum initial perturbation", "", nlsP->imin, &nlsP->imin,NULL);CHKERRQ(ierr); 967 ierr = PetscOptionsReal("-tao_nls_imax", "maximum initial perturbation", "", nlsP->imax, &nlsP->imax,NULL);CHKERRQ(ierr); 968 ierr = PetscOptionsReal("-tao_nls_imfac", "initial merit factor", "", nlsP->imfac, &nlsP->imfac,NULL);CHKERRQ(ierr); 969 ierr = PetscOptionsReal("-tao_nls_pmin", "minimum perturbation", "", nlsP->pmin, &nlsP->pmin,NULL);CHKERRQ(ierr); 970 ierr = PetscOptionsReal("-tao_nls_pmax", "maximum perturbation", "", nlsP->pmax, &nlsP->pmax,NULL);CHKERRQ(ierr); 971 ierr = PetscOptionsReal("-tao_nls_pgfac", "growth factor", "", nlsP->pgfac, &nlsP->pgfac,NULL);CHKERRQ(ierr); 972 ierr = PetscOptionsReal("-tao_nls_psfac", "shrink factor", "", nlsP->psfac, &nlsP->psfac,NULL);CHKERRQ(ierr); 973 ierr = PetscOptionsReal("-tao_nls_pmgfac", "merit growth factor", "", nlsP->pmgfac, &nlsP->pmgfac,NULL);CHKERRQ(ierr); 974 ierr = PetscOptionsReal("-tao_nls_pmsfac", "merit shrink factor", "", nlsP->pmsfac, &nlsP->pmsfac,NULL);CHKERRQ(ierr); 975 ierr = PetscOptionsReal("-tao_nls_eta1", "poor steplength; reduce radius", "", nlsP->eta1, &nlsP->eta1,NULL);CHKERRQ(ierr); 976 ierr = PetscOptionsReal("-tao_nls_eta2", "reasonable steplength; leave radius alone", "", nlsP->eta2, &nlsP->eta2,NULL);CHKERRQ(ierr); 977 ierr = PetscOptionsReal("-tao_nls_eta3", "good steplength; increase radius", "", nlsP->eta3, &nlsP->eta3,NULL);CHKERRQ(ierr); 978 ierr = PetscOptionsReal("-tao_nls_eta4", "excellent steplength; greatly increase radius", "", nlsP->eta4, &nlsP->eta4,NULL);CHKERRQ(ierr); 979 ierr = PetscOptionsReal("-tao_nls_alpha1", "", "", nlsP->alpha1, &nlsP->alpha1,NULL);CHKERRQ(ierr); 980 ierr = PetscOptionsReal("-tao_nls_alpha2", "", "", nlsP->alpha2, &nlsP->alpha2,NULL);CHKERRQ(ierr); 981 ierr = PetscOptionsReal("-tao_nls_alpha3", "", "", nlsP->alpha3, &nlsP->alpha3,NULL);CHKERRQ(ierr); 982 ierr = PetscOptionsReal("-tao_nls_alpha4", "", "", nlsP->alpha4, &nlsP->alpha4,NULL);CHKERRQ(ierr); 983 ierr = PetscOptionsReal("-tao_nls_alpha5", "", "", nlsP->alpha5, &nlsP->alpha5,NULL);CHKERRQ(ierr); 984 ierr = PetscOptionsReal("-tao_nls_nu1", "poor steplength; reduce radius", "", nlsP->nu1, &nlsP->nu1,NULL);CHKERRQ(ierr); 985 ierr = PetscOptionsReal("-tao_nls_nu2", "reasonable steplength; leave radius alone", "", nlsP->nu2, &nlsP->nu2,NULL);CHKERRQ(ierr); 986 ierr = PetscOptionsReal("-tao_nls_nu3", "good steplength; increase radius", "", nlsP->nu3, &nlsP->nu3,NULL);CHKERRQ(ierr); 987 ierr = PetscOptionsReal("-tao_nls_nu4", "excellent steplength; greatly increase radius", "", nlsP->nu4, &nlsP->nu4,NULL);CHKERRQ(ierr); 988 ierr = PetscOptionsReal("-tao_nls_omega1", "", "", nlsP->omega1, &nlsP->omega1,NULL);CHKERRQ(ierr); 989 ierr = PetscOptionsReal("-tao_nls_omega2", "", "", nlsP->omega2, &nlsP->omega2,NULL);CHKERRQ(ierr); 990 ierr = PetscOptionsReal("-tao_nls_omega3", "", "", nlsP->omega3, &nlsP->omega3,NULL);CHKERRQ(ierr); 991 ierr = PetscOptionsReal("-tao_nls_omega4", "", "", nlsP->omega4, &nlsP->omega4,NULL);CHKERRQ(ierr); 992 ierr = PetscOptionsReal("-tao_nls_omega5", "", "", nlsP->omega5, &nlsP->omega5,NULL);CHKERRQ(ierr); 993 ierr = PetscOptionsReal("-tao_nls_mu1_i", "", "", nlsP->mu1_i, &nlsP->mu1_i,NULL);CHKERRQ(ierr); 994 ierr = PetscOptionsReal("-tao_nls_mu2_i", "", "", nlsP->mu2_i, &nlsP->mu2_i,NULL);CHKERRQ(ierr); 995 ierr = PetscOptionsReal("-tao_nls_gamma1_i", "", "", nlsP->gamma1_i, &nlsP->gamma1_i,NULL);CHKERRQ(ierr); 996 ierr = PetscOptionsReal("-tao_nls_gamma2_i", "", "", nlsP->gamma2_i, &nlsP->gamma2_i,NULL);CHKERRQ(ierr); 997 ierr = PetscOptionsReal("-tao_nls_gamma3_i", "", "", nlsP->gamma3_i, &nlsP->gamma3_i,NULL);CHKERRQ(ierr); 998 ierr = PetscOptionsReal("-tao_nls_gamma4_i", "", "", nlsP->gamma4_i, &nlsP->gamma4_i,NULL);CHKERRQ(ierr); 999 ierr = PetscOptionsReal("-tao_nls_theta_i", "", "", nlsP->theta_i, &nlsP->theta_i,NULL);CHKERRQ(ierr); 1000 ierr = PetscOptionsReal("-tao_nls_mu1", "", "", nlsP->mu1, &nlsP->mu1,NULL);CHKERRQ(ierr); 1001 ierr = PetscOptionsReal("-tao_nls_mu2", "", "", nlsP->mu2, &nlsP->mu2,NULL);CHKERRQ(ierr); 1002 ierr = PetscOptionsReal("-tao_nls_gamma1", "", "", nlsP->gamma1, &nlsP->gamma1,NULL);CHKERRQ(ierr); 1003 ierr = PetscOptionsReal("-tao_nls_gamma2", "", "", nlsP->gamma2, &nlsP->gamma2,NULL);CHKERRQ(ierr); 1004 ierr = PetscOptionsReal("-tao_nls_gamma3", "", "", nlsP->gamma3, &nlsP->gamma3,NULL);CHKERRQ(ierr); 1005 ierr = PetscOptionsReal("-tao_nls_gamma4", "", "", nlsP->gamma4, &nlsP->gamma4,NULL);CHKERRQ(ierr); 1006 ierr = PetscOptionsReal("-tao_nls_theta", "", "", nlsP->theta, &nlsP->theta,NULL);CHKERRQ(ierr); 1007 ierr = PetscOptionsReal("-tao_nls_min_radius", "lower bound on initial radius", "", nlsP->min_radius, &nlsP->min_radius,NULL);CHKERRQ(ierr); 1008 ierr = PetscOptionsReal("-tao_nls_max_radius", "upper bound on radius", "", nlsP->max_radius, &nlsP->max_radius,NULL);CHKERRQ(ierr); 1009 ierr = PetscOptionsReal("-tao_nls_epsilon", "tolerance used when computing actual and predicted reduction", "", nlsP->epsilon, &nlsP->epsilon,NULL);CHKERRQ(ierr); 1010 ierr = PetscOptionsTail();CHKERRQ(ierr); 1011 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1012 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 1013 PetscFunctionReturn(0); 1014 } 1015 1016 1017 /*------------------------------------------------------------*/ 1018 #undef __FUNCT__ 1019 #define __FUNCT__ "TaoView_NLS" 1020 static PetscErrorCode TaoView_NLS(Tao tao, PetscViewer viewer) 1021 { 1022 TAO_NLS *nlsP = (TAO_NLS *)tao->data; 1023 PetscInt nrejects; 1024 PetscBool isascii; 1025 PetscErrorCode ierr; 1026 1027 PetscFunctionBegin; 1028 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1029 if (isascii) { 1030 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1031 if (NLS_PC_BFGS == nlsP->pc_type && nlsP->M) { 1032 ierr = MatLMVMGetRejects(nlsP->M,&nrejects);CHKERRQ(ierr); 1033 ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr); 1034 } 1035 ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", nlsP->newt);CHKERRQ(ierr); 1036 ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", nlsP->bfgs);CHKERRQ(ierr); 1037 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", nlsP->sgrad);CHKERRQ(ierr); 1038 ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", nlsP->grad);CHKERRQ(ierr); 1039 1040 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp atol: %D\n", nlsP->ksp_atol);CHKERRQ(ierr); 1041 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp rtol: %D\n", nlsP->ksp_rtol);CHKERRQ(ierr); 1042 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp ctol: %D\n", nlsP->ksp_ctol);CHKERRQ(ierr); 1043 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp negc: %D\n", nlsP->ksp_negc);CHKERRQ(ierr); 1044 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp dtol: %D\n", nlsP->ksp_dtol);CHKERRQ(ierr); 1045 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp iter: %D\n", nlsP->ksp_iter);CHKERRQ(ierr); 1046 ierr = PetscViewerASCIIPrintf(viewer, "nls ksp othr: %D\n", nlsP->ksp_othr);CHKERRQ(ierr); 1047 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1048 } 1049 PetscFunctionReturn(0); 1050 } 1051 1052 /* ---------------------------------------------------------- */ 1053 /*MC 1054 TAONLS - Newton's method with linesearch for unconstrained minimization. 1055 At each iteration, the Newton line search method solves the symmetric 1056 system of equations to obtain the step diretion dk: 1057 Hk dk = -gk 1058 a More-Thuente line search is applied on the direction dk to approximately 1059 solve 1060 min_t f(xk + t d_k) 1061 1062 Options Database Keys: 1063 + -tao_nls_ksp_type - "cg","nash","stcg","gltr","petsc" 1064 . -tao_nls_pc_type - "none","ahess","bfgs","petsc" 1065 . -tao_nls_bfgs_scale_type - "ahess","phess","bfgs" 1066 . -tao_nls_init_type - "constant","direction","interpolation" 1067 . -tao_nls_update_type - "step","direction","interpolation" 1068 . -tao_nls_sval - perturbation starting value 1069 . -tao_nls_imin - minimum initial perturbation 1070 . -tao_nls_imax - maximum initial perturbation 1071 . -tao_nls_pmin - minimum perturbation 1072 . -tao_nls_pmax - maximum perturbation 1073 . -tao_nls_pgfac - growth factor 1074 . -tao_nls_psfac - shrink factor 1075 . -tao_nls_imfac - initial merit factor 1076 . -tao_nls_pmgfac - merit growth factor 1077 . -tao_nls_pmsfac - merit shrink factor 1078 . -tao_nls_eta1 - poor steplength; reduce radius 1079 . -tao_nls_eta2 - reasonable steplength; leave radius 1080 . -tao_nls_eta3 - good steplength; increase readius 1081 . -tao_nls_eta4 - excellent steplength; greatly increase radius 1082 . -tao_nls_alpha1 - alpha1 reduction 1083 . -tao_nls_alpha2 - alpha2 reduction 1084 . -tao_nls_alpha3 - alpha3 reduction 1085 . -tao_nls_alpha4 - alpha4 reduction 1086 . -tao_nls_alpha - alpha5 reduction 1087 . -tao_nls_mu1 - mu1 interpolation update 1088 . -tao_nls_mu2 - mu2 interpolation update 1089 . -tao_nls_gamma1 - gamma1 interpolation update 1090 . -tao_nls_gamma2 - gamma2 interpolation update 1091 . -tao_nls_gamma3 - gamma3 interpolation update 1092 . -tao_nls_gamma4 - gamma4 interpolation update 1093 . -tao_nls_theta - theta interpolation update 1094 . -tao_nls_omega1 - omega1 step update 1095 . -tao_nls_omega2 - omega2 step update 1096 . -tao_nls_omega3 - omega3 step update 1097 . -tao_nls_omega4 - omega4 step update 1098 . -tao_nls_omega5 - omega5 step update 1099 . -tao_nls_mu1_i - mu1 interpolation init factor 1100 . -tao_nls_mu2_i - mu2 interpolation init factor 1101 . -tao_nls_gamma1_i - gamma1 interpolation init factor 1102 . -tao_nls_gamma2_i - gamma2 interpolation init factor 1103 . -tao_nls_gamma3_i - gamma3 interpolation init factor 1104 . -tao_nls_gamma4_i - gamma4 interpolation init factor 1105 - -tao_nls_theta_i - theta interpolation init factor 1106 1107 Level: beginner 1108 M*/ 1109 1110 #undef __FUNCT__ 1111 #define __FUNCT__ "TaoCreate_NLS" 1112 PETSC_EXTERN PetscErrorCode TaoCreate_NLS(Tao tao) 1113 { 1114 TAO_NLS *nlsP; 1115 const char *morethuente_type = TAOLINESEARCHMT; 1116 PetscErrorCode ierr; 1117 1118 PetscFunctionBegin; 1119 ierr = PetscNewLog(tao,&nlsP);CHKERRQ(ierr); 1120 1121 tao->ops->setup = TaoSetUp_NLS; 1122 tao->ops->solve = TaoSolve_NLS; 1123 tao->ops->view = TaoView_NLS; 1124 tao->ops->setfromoptions = TaoSetFromOptions_NLS; 1125 tao->ops->destroy = TaoDestroy_NLS; 1126 1127 tao->max_it = 50; 1128 #if defined(PETSC_USE_REAL_SINGLE) 1129 tao->fatol = 1e-5; 1130 tao->frtol = 1e-5; 1131 #else 1132 tao->fatol = 1e-10; 1133 tao->frtol = 1e-10; 1134 #endif 1135 tao->data = (void*)nlsP; 1136 tao->trust0 = 100.0; 1137 1138 nlsP->sval = 0.0; 1139 nlsP->imin = 1.0e-4; 1140 nlsP->imax = 1.0e+2; 1141 nlsP->imfac = 1.0e-1; 1142 1143 nlsP->pmin = 1.0e-12; 1144 nlsP->pmax = 1.0e+2; 1145 nlsP->pgfac = 1.0e+1; 1146 nlsP->psfac = 4.0e-1; 1147 nlsP->pmgfac = 1.0e-1; 1148 nlsP->pmsfac = 1.0e-1; 1149 1150 /* Default values for trust-region radius update based on steplength */ 1151 nlsP->nu1 = 0.25; 1152 nlsP->nu2 = 0.50; 1153 nlsP->nu3 = 1.00; 1154 nlsP->nu4 = 1.25; 1155 1156 nlsP->omega1 = 0.25; 1157 nlsP->omega2 = 0.50; 1158 nlsP->omega3 = 1.00; 1159 nlsP->omega4 = 2.00; 1160 nlsP->omega5 = 4.00; 1161 1162 /* Default values for trust-region radius update based on reduction */ 1163 nlsP->eta1 = 1.0e-4; 1164 nlsP->eta2 = 0.25; 1165 nlsP->eta3 = 0.50; 1166 nlsP->eta4 = 0.90; 1167 1168 nlsP->alpha1 = 0.25; 1169 nlsP->alpha2 = 0.50; 1170 nlsP->alpha3 = 1.00; 1171 nlsP->alpha4 = 2.00; 1172 nlsP->alpha5 = 4.00; 1173 1174 /* Default values for trust-region radius update based on interpolation */ 1175 nlsP->mu1 = 0.10; 1176 nlsP->mu2 = 0.50; 1177 1178 nlsP->gamma1 = 0.25; 1179 nlsP->gamma2 = 0.50; 1180 nlsP->gamma3 = 2.00; 1181 nlsP->gamma4 = 4.00; 1182 1183 nlsP->theta = 0.05; 1184 1185 /* Default values for trust region initialization based on interpolation */ 1186 nlsP->mu1_i = 0.35; 1187 nlsP->mu2_i = 0.50; 1188 1189 nlsP->gamma1_i = 0.0625; 1190 nlsP->gamma2_i = 0.5; 1191 nlsP->gamma3_i = 2.0; 1192 nlsP->gamma4_i = 5.0; 1193 1194 nlsP->theta_i = 0.25; 1195 1196 /* Remaining parameters */ 1197 nlsP->min_radius = 1.0e-10; 1198 nlsP->max_radius = 1.0e10; 1199 nlsP->epsilon = 1.0e-6; 1200 1201 nlsP->ksp_type = NLS_KSP_STCG; 1202 nlsP->pc_type = NLS_PC_BFGS; 1203 nlsP->bfgs_scale_type = BFGS_SCALE_PHESS; 1204 nlsP->init_type = NLS_INIT_INTERPOLATION; 1205 nlsP->update_type = NLS_UPDATE_STEP; 1206 1207 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr); 1208 ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 1209 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 1210 1211 /* Set linear solver to default for symmetric matrices */ 1212 ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr); 1213 PetscFunctionReturn(0); 1214 } 1215 1216 #undef __FUNCT__ 1217 #define __FUNCT__ "MatLMVMSolveShell" 1218 static PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x) 1219 { 1220 PetscErrorCode ierr; 1221 Mat M; 1222 1223 PetscFunctionBegin; 1224 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1225 PetscValidHeaderSpecific(b,VEC_CLASSID,2); 1226 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 1227 ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr); 1228 ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr); 1229 PetscFunctionReturn(0); 1230 } 1231