1 #define PETSCSNES_DLL 2 3 #include "../src/snes/impls/ls/ls.h" 4 5 /* 6 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 7 || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that 8 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 9 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 10 */ 11 #undef __FUNCT__ 12 #define __FUNCT__ "SNESLSCheckLocalMin_Private" 13 PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 14 { 15 PetscReal a1; 16 PetscErrorCode ierr; 17 PetscTruth hastranspose; 18 19 PetscFunctionBegin; 20 *ismin = PETSC_FALSE; 21 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 22 if (hastranspose) { 23 /* Compute || J^T F|| */ 24 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 25 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 26 ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 27 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 28 } else { 29 Vec work; 30 PetscScalar result; 31 PetscReal wnorm; 32 33 ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 34 ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 35 ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 36 ierr = MatMult(A,W,work);CHKERRQ(ierr); 37 ierr = VecDot(F,work,&result);CHKERRQ(ierr); 38 ierr = VecDestroy(work);CHKERRQ(ierr); 39 a1 = PetscAbsScalar(result)/(fnorm*wnorm); 40 ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 41 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 42 } 43 PetscFunctionReturn(0); 44 } 45 46 /* 47 Checks if J^T(F - J*X) = 0 48 */ 49 #undef __FUNCT__ 50 #define __FUNCT__ "SNESLSCheckResidual_Private" 51 PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 52 { 53 PetscReal a1,a2; 54 PetscErrorCode ierr; 55 PetscTruth hastranspose; 56 57 PetscFunctionBegin; 58 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 59 if (hastranspose) { 60 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 61 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 62 63 /* Compute || J^T W|| */ 64 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 65 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 66 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 67 if (a1) { 68 ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 69 } 70 } 71 PetscFunctionReturn(0); 72 } 73 74 /* -------------------------------------------------------------------- 75 76 This file implements a truncated Newton method with a line search, 77 for solving a system of nonlinear equations, using the KSP, Vec, 78 and Mat interfaces for linear solvers, vectors, and matrices, 79 respectively. 80 81 The following basic routines are required for each nonlinear solver: 82 SNESCreate_XXX() - Creates a nonlinear solver context 83 SNESSetFromOptions_XXX() - Sets runtime options 84 SNESSolve_XXX() - Solves the nonlinear system 85 SNESDestroy_XXX() - Destroys the nonlinear solver context 86 The suffix "_XXX" denotes a particular implementation, in this case 87 we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 88 systems of nonlinear equations with a line search (LS) method. 89 These routines are actually called via the common user interface 90 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 91 SNESDestroy(), so the application code interface remains identical 92 for all nonlinear solvers. 93 94 Another key routine is: 95 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 96 by setting data structures and options. The interface routine SNESSetUp() 97 is not usually called directly by the user, but instead is called by 98 SNESSolve() if necessary. 99 100 Additional basic routines are: 101 SNESView_XXX() - Prints details of runtime options that 102 have actually been used. 103 These are called by application codes via the interface routines 104 SNESView(). 105 106 The various types of solvers (preconditioners, Krylov subspace methods, 107 nonlinear solvers, timesteppers) are all organized similarly, so the 108 above description applies to these categories also. 109 110 -------------------------------------------------------------------- */ 111 /* 112 SNESSolve_LS - Solves a nonlinear system with a truncated Newton 113 method with a line search. 114 115 Input Parameters: 116 . snes - the SNES context 117 118 Output Parameter: 119 . outits - number of iterations until termination 120 121 Application Interface Routine: SNESSolve() 122 123 Notes: 124 This implements essentially a truncated Newton method with a 125 line search. By default a cubic backtracking line search 126 is employed, as described in the text "Numerical Methods for 127 Unconstrained Optimization and Nonlinear Equations" by Dennis 128 and Schnabel. 129 */ 130 #undef __FUNCT__ 131 #define __FUNCT__ "SNESSolve_LS" 132 PetscErrorCode SNESSolve_LS(SNES snes) 133 { 134 SNES_LS *neP = (SNES_LS*)snes->data; 135 PetscErrorCode ierr; 136 PetscInt maxits,i,lits; 137 PetscTruth lssucceed; 138 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 139 PetscReal fnorm,gnorm,xnorm=0,ynorm; 140 Vec Y,X,F,G,W; 141 KSPConvergedReason kspreason; 142 143 PetscFunctionBegin; 144 snes->numFailures = 0; 145 snes->numLinearSolveFailures = 0; 146 snes->reason = SNES_CONVERGED_ITERATING; 147 148 maxits = snes->max_its; /* maximum number of iterations */ 149 X = snes->vec_sol; /* solution vector */ 150 F = snes->vec_func; /* residual vector */ 151 Y = snes->work[0]; /* work vectors */ 152 G = snes->work[1]; 153 W = snes->work[2]; 154 155 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 156 snes->iter = 0; 157 snes->norm = 0; 158 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 159 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 160 if (snes->domainerror) { 161 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 162 PetscFunctionReturn(0); 163 } 164 ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 165 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 166 ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 167 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 168 if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 169 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 170 snes->norm = fnorm; 171 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 172 SNESLogConvHistory(snes,fnorm,0); 173 SNESMonitor(snes,0,fnorm); 174 175 /* set parameter for default relative tolerance convergence test */ 176 snes->ttol = fnorm*snes->rtol; 177 /* test convergence */ 178 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 179 if (snes->reason) PetscFunctionReturn(0); 180 181 for (i=0; i<maxits; i++) { 182 183 /* Call general purpose update function */ 184 if (snes->ops->update) { 185 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 186 } 187 188 /* Solve J Y = F, where J is Jacobian matrix */ 189 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 190 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 191 ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 192 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 193 if (kspreason < 0) { 194 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 195 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 196 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 197 break; 198 } 199 } 200 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 201 snes->linear_its += lits; 202 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 203 204 if (neP->precheckstep) { 205 PetscTruth changed_y = PETSC_FALSE; 206 ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 207 } 208 209 if (PetscLogPrintInfo){ 210 ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 211 } 212 213 /* Compute a (scaled) negative update in the line search routine: 214 Y <- X - lambda*Y 215 and evaluate G = function(Y) (depends on the line search). 216 */ 217 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 218 ynorm = 1; gnorm = fnorm; 219 ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 220 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 221 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 222 if (snes->domainerror) { 223 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 224 PetscFunctionReturn(0); 225 } 226 if (!lssucceed) { 227 if (++snes->numFailures >= snes->maxFailures) { 228 PetscTruth ismin; 229 snes->reason = SNES_DIVERGED_LS_FAILURE; 230 ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 231 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 232 break; 233 } 234 } 235 /* Update function and solution vectors */ 236 fnorm = gnorm; 237 ierr = VecCopy(G,F);CHKERRQ(ierr); 238 ierr = VecCopy(W,X);CHKERRQ(ierr); 239 /* Monitor convergence */ 240 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 241 snes->iter = i+1; 242 snes->norm = fnorm; 243 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 244 SNESLogConvHistory(snes,snes->norm,lits); 245 SNESMonitor(snes,snes->iter,snes->norm); 246 /* Test for convergence, xnorm = || X || */ 247 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 248 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 249 if (snes->reason) break; 250 } 251 if (i == maxits) { 252 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 253 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 254 } 255 PetscFunctionReturn(0); 256 } 257 /* -------------------------------------------------------------------------- */ 258 /* 259 SNESSetUp_LS - Sets up the internal data structures for the later use 260 of the SNESLS nonlinear solver. 261 262 Input Parameter: 263 . snes - the SNES context 264 . x - the solution vector 265 266 Application Interface Routine: SNESSetUp() 267 268 Notes: 269 For basic use of the SNES solvers, the user need not explicitly call 270 SNESSetUp(), since these actions will automatically occur during 271 the call to SNESSolve(). 272 */ 273 #undef __FUNCT__ 274 #define __FUNCT__ "SNESSetUp_LS" 275 PetscErrorCode SNESSetUp_LS(SNES snes) 276 { 277 PetscErrorCode ierr; 278 279 PetscFunctionBegin; 280 if (!snes->vec_sol_update) { 281 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 282 ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 283 } 284 if (!snes->work) { 285 snes->nwork = 3; 286 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 287 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 288 } 289 PetscFunctionReturn(0); 290 } 291 /* -------------------------------------------------------------------------- */ 292 /* 293 SNESDestroy_LS - Destroys the private SNES_LS context that was created 294 with SNESCreate_LS(). 295 296 Input Parameter: 297 . snes - the SNES context 298 299 Application Interface Routine: SNESDestroy() 300 */ 301 #undef __FUNCT__ 302 #define __FUNCT__ "SNESDestroy_LS" 303 PetscErrorCode SNESDestroy_LS(SNES snes) 304 { 305 PetscErrorCode ierr; 306 307 PetscFunctionBegin; 308 if (snes->vec_sol_update) { 309 ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr); 310 snes->vec_sol_update = PETSC_NULL; 311 } 312 if (snes->nwork) { 313 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 314 snes->nwork = 0; 315 snes->work = PETSC_NULL; 316 } 317 ierr = PetscFree(snes->data);CHKERRQ(ierr); 318 319 /* clear composed functions */ 320 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 321 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 322 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 323 324 PetscFunctionReturn(0); 325 } 326 /* -------------------------------------------------------------------------- */ 327 #undef __FUNCT__ 328 #define __FUNCT__ "SNESLineSearchNo" 329 330 /*@C 331 SNESLineSearchNo - This routine is not a line search at all; 332 it simply uses the full Newton step. Thus, this routine is intended 333 to serve as a template and is not recommended for general use. 334 335 Collective on SNES and Vec 336 337 Input Parameters: 338 + snes - nonlinear context 339 . lsctx - optional context for line search (not used here) 340 . x - current iterate 341 . f - residual evaluated at x 342 . y - search direction 343 . w - work vector 344 . fnorm - 2-norm of f 345 - xnorm - norm of x if known, otherwise 0 346 347 Output Parameters: 348 + g - residual evaluated at new iterate y 349 . w - new iterate 350 . gnorm - 2-norm of g 351 . ynorm - 2-norm of search length 352 - flag - PETSC_TRUE on success, PETSC_FALSE on failure 353 354 Options Database Key: 355 . -snes_ls basic - Activates SNESLineSearchNo() 356 357 Level: advanced 358 359 .keywords: SNES, nonlinear, line search, cubic 360 361 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 362 SNESLineSearchSet(), SNESLineSearchNoNorms() 363 @*/ 364 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 365 { 366 PetscErrorCode ierr; 367 SNES_LS *neP = (SNES_LS*)snes->data; 368 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 369 370 PetscFunctionBegin; 371 *flag = PETSC_TRUE; 372 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 373 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 374 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 375 if (neP->postcheckstep) { 376 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 377 } 378 if (changed_y) { 379 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 380 } 381 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 382 if (!snes->domainerror) { 383 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 384 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 385 } 386 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 387 PetscFunctionReturn(0); 388 } 389 /* -------------------------------------------------------------------------- */ 390 #undef __FUNCT__ 391 #define __FUNCT__ "SNESLineSearchNoNorms" 392 393 /*@C 394 SNESLineSearchNoNorms - This routine is not a line search at 395 all; it simply uses the full Newton step. This version does not 396 even compute the norm of the function or search direction; this 397 is intended only when you know the full step is fine and are 398 not checking for convergence of the nonlinear iteration (for 399 example, you are running always for a fixed number of Newton steps). 400 401 Collective on SNES and Vec 402 403 Input Parameters: 404 + snes - nonlinear context 405 . lsctx - optional context for line search (not used here) 406 . x - current iterate 407 . f - residual evaluated at x 408 . y - search direction 409 . w - work vector 410 . fnorm - 2-norm of f 411 - xnorm - norm of x if known, otherwise 0 412 413 Output Parameters: 414 + g - residual evaluated at new iterate y 415 . w - new iterate 416 . gnorm - not changed 417 . ynorm - not changed 418 - flag - set to PETSC_TRUE indicating a successful line search 419 420 Options Database Key: 421 . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 422 423 Notes: 424 SNESLineSearchNoNorms() must be used in conjunction with 425 either the options 426 $ -snes_no_convergence_test -snes_max_it <its> 427 or alternatively a user-defined custom test set via 428 SNESSetConvergenceTest(); or a -snes_max_it of 1, 429 otherwise, the SNES solver will generate an error. 430 431 During the final iteration this will not evaluate the function at 432 the solution point. This is to save a function evaluation while 433 using pseudo-timestepping. 434 435 The residual norms printed by monitoring routines such as 436 SNESMonitorDefault() (as activated via -snes_monitor) will not be 437 correct, since they are not computed. 438 439 Level: advanced 440 441 .keywords: SNES, nonlinear, line search, cubic 442 443 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 444 SNESLineSearchSet(), SNESLineSearchNo() 445 @*/ 446 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 447 { 448 PetscErrorCode ierr; 449 SNES_LS *neP = (SNES_LS*)snes->data; 450 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 451 452 PetscFunctionBegin; 453 *flag = PETSC_TRUE; 454 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 455 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 456 if (neP->postcheckstep) { 457 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 458 } 459 if (changed_y) { 460 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 461 } 462 463 /* don't evaluate function the last time through */ 464 if (snes->iter < snes->max_its-1) { 465 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 466 } 467 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 468 PetscFunctionReturn(0); 469 } 470 /* -------------------------------------------------------------------------- */ 471 #undef __FUNCT__ 472 #define __FUNCT__ "SNESLineSearchCubic" 473 /*@C 474 SNESLineSearchCubic - Performs a cubic line search (default line search method). 475 476 Collective on SNES 477 478 Input Parameters: 479 + snes - nonlinear context 480 . lsctx - optional context for line search (not used here) 481 . x - current iterate 482 . f - residual evaluated at x 483 . y - search direction 484 . w - work vector 485 . fnorm - 2-norm of f 486 - xnorm - norm of x if known, otherwise 0 487 488 Output Parameters: 489 + g - residual evaluated at new iterate y 490 . w - new iterate 491 . gnorm - 2-norm of g 492 . ynorm - 2-norm of search length 493 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 494 495 Options Database Key: 496 + -snes_ls cubic - Activates SNESLineSearchCubic() 497 . -snes_ls_alpha <alpha> - Sets alpha 498 . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 499 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 500 501 502 Notes: 503 This line search is taken from "Numerical Methods for Unconstrained 504 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 505 506 Level: advanced 507 508 .keywords: SNES, nonlinear, line search, cubic 509 510 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 511 @*/ 512 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 513 { 514 /* 515 Note that for line search purposes we work with with the related 516 minimization problem: 517 min z(x): R^n -> R, 518 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 519 */ 520 521 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 522 PetscReal minlambda,lambda,lambdatemp; 523 #if defined(PETSC_USE_COMPLEX) 524 PetscScalar cinitslope; 525 #endif 526 PetscErrorCode ierr; 527 PetscInt count; 528 SNES_LS *neP = (SNES_LS*)snes->data; 529 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 530 531 PetscFunctionBegin; 532 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 533 *flag = PETSC_TRUE; 534 535 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 536 if (!*ynorm) { 537 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 538 *gnorm = fnorm; 539 ierr = VecCopy(x,w);CHKERRQ(ierr); 540 ierr = VecCopy(f,g);CHKERRQ(ierr); 541 *flag = PETSC_FALSE; 542 goto theend1; 543 } 544 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 545 ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 546 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 547 *ynorm = neP->maxstep; 548 } 549 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 550 minlambda = neP->minlambda/rellength; 551 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 552 #if defined(PETSC_USE_COMPLEX) 553 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 554 initslope = PetscRealPart(cinitslope); 555 #else 556 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 557 #endif 558 if (initslope > 0.0) initslope = -initslope; 559 if (initslope == 0.0) initslope = -1.0; 560 561 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 562 if (snes->nfuncs >= snes->max_funcs) { 563 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 564 *flag = PETSC_FALSE; 565 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 566 goto theend1; 567 } 568 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 569 if (snes->domainerror) { 570 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 571 PetscFunctionReturn(0); 572 } 573 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 574 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 575 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 576 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 577 ierr = PetscInfo2(snes,"Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 578 goto theend1; 579 } 580 581 /* Fit points with quadratic */ 582 lambda = 1.0; 583 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 584 lambdaprev = lambda; 585 gnormprev = *gnorm; 586 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 587 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 588 else lambda = lambdatemp; 589 590 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 591 if (snes->nfuncs >= snes->max_funcs) { 592 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 593 *flag = PETSC_FALSE; 594 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 595 goto theend1; 596 } 597 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 598 if (snes->domainerror) { 599 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 600 PetscFunctionReturn(0); 601 } 602 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 603 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 604 ierr = PetscInfo1(snes,"gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 605 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 606 ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 607 goto theend1; 608 } 609 610 /* Fit points with cubic */ 611 count = 1; 612 while (PETSC_TRUE) { 613 if (lambda <= minlambda) { 614 ierr = PetscInfo1(snes,"Unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 615 ierr = PetscInfo6(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr); 616 *flag = PETSC_FALSE; 617 break; 618 } 619 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 620 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 621 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 622 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 623 d = b*b - 3*a*initslope; 624 if (d < 0.0) d = 0.0; 625 if (a == 0.0) { 626 lambdatemp = -initslope/(2.0*b); 627 } else { 628 lambdatemp = (-b + sqrt(d))/(3.0*a); 629 } 630 lambdaprev = lambda; 631 gnormprev = *gnorm; 632 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 633 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 634 else lambda = lambdatemp; 635 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 636 if (snes->nfuncs >= snes->max_funcs) { 637 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 638 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 639 *flag = PETSC_FALSE; 640 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 641 break; 642 } 643 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 644 if (snes->domainerror) { 645 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 646 PetscFunctionReturn(0); 647 } 648 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 649 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 650 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 651 ierr = PetscInfo2(snes,"Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 652 break; 653 } else { 654 ierr = PetscInfo2(snes,"Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 655 } 656 count++; 657 } 658 theend1: 659 /* Optional user-defined check for line search step validity */ 660 if (neP->postcheckstep && *flag) { 661 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 662 if (changed_y) { 663 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 664 } 665 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 666 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 667 if (snes->domainerror) { 668 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 669 PetscFunctionReturn(0); 670 } 671 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 672 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 673 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 674 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 675 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 676 } 677 } 678 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 679 PetscFunctionReturn(0); 680 } 681 /* -------------------------------------------------------------------------- */ 682 #undef __FUNCT__ 683 #define __FUNCT__ "SNESLineSearchQuadratic" 684 /*@C 685 SNESLineSearchQuadratic - Performs a quadratic line search. 686 687 Collective on SNES and Vec 688 689 Input Parameters: 690 + snes - the SNES context 691 . lsctx - optional context for line search (not used here) 692 . x - current iterate 693 . f - residual evaluated at x 694 . y - search direction 695 . w - work vector 696 . fnorm - 2-norm of f 697 - xnorm - norm of x if known, otherwise 0 698 699 Output Parameters: 700 + g - residual evaluated at new iterate w 701 . w - new iterate (x + lambda*y) 702 . gnorm - 2-norm of g 703 . ynorm - 2-norm of search length 704 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 705 706 Options Database Keys: 707 + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 708 . -snes_ls_alpha <alpha> - Sets alpha 709 . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 710 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 711 712 Notes: 713 Use SNESLineSearchSet() to set this routine within the SNESLS method. 714 715 Level: advanced 716 717 .keywords: SNES, nonlinear, quadratic, line search 718 719 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 720 @*/ 721 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 722 { 723 /* 724 Note that for line search purposes we work with with the related 725 minimization problem: 726 min z(x): R^n -> R, 727 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 728 */ 729 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 730 #if defined(PETSC_USE_COMPLEX) 731 PetscScalar cinitslope; 732 #endif 733 PetscErrorCode ierr; 734 PetscInt count; 735 SNES_LS *neP = (SNES_LS*)snes->data; 736 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 737 738 PetscFunctionBegin; 739 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 740 *flag = PETSC_TRUE; 741 742 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 743 if (*ynorm == 0.0) { 744 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 745 *gnorm = fnorm; 746 ierr = VecCopy(x,w);CHKERRQ(ierr); 747 ierr = VecCopy(f,g);CHKERRQ(ierr); 748 *flag = PETSC_FALSE; 749 goto theend2; 750 } 751 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 752 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 753 *ynorm = neP->maxstep; 754 } 755 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 756 minlambda = neP->minlambda/rellength; 757 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 758 #if defined(PETSC_USE_COMPLEX) 759 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 760 initslope = PetscRealPart(cinitslope); 761 #else 762 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 763 #endif 764 if (initslope > 0.0) initslope = -initslope; 765 if (initslope == 0.0) initslope = -1.0; 766 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 767 768 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 769 if (snes->nfuncs >= snes->max_funcs) { 770 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 771 *flag = PETSC_FALSE; 772 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 773 goto theend2; 774 } 775 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 776 if (snes->domainerror) { 777 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 778 PetscFunctionReturn(0); 779 } 780 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 781 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 782 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 783 ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 784 goto theend2; 785 } 786 787 /* Fit points with quadratic */ 788 lambda = 1.0; 789 count = 1; 790 while (PETSC_TRUE) { 791 if (lambda <= minlambda) { /* bad luck; use full step */ 792 ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 793 ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 794 ierr = VecCopy(x,w);CHKERRQ(ierr); 795 *flag = PETSC_FALSE; 796 break; 797 } 798 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 799 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 800 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 801 else lambda = lambdatemp; 802 803 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 804 if (snes->nfuncs >= snes->max_funcs) { 805 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 806 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 807 *flag = PETSC_FALSE; 808 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 809 break; 810 } 811 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 812 if (snes->domainerror) { 813 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 814 PetscFunctionReturn(0); 815 } 816 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 817 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 818 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 819 ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 820 break; 821 } 822 count++; 823 } 824 theend2: 825 /* Optional user-defined check for line search step validity */ 826 if (neP->postcheckstep) { 827 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 828 if (changed_y) { 829 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 830 } 831 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 832 ierr = SNESComputeFunction(snes,w,g); 833 if (snes->domainerror) { 834 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 835 PetscFunctionReturn(0); 836 } 837 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 838 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 839 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 840 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 841 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 842 } 843 } 844 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 845 PetscFunctionReturn(0); 846 } 847 848 /* -------------------------------------------------------------------------- */ 849 #undef __FUNCT__ 850 #define __FUNCT__ "SNESLineSearchSet" 851 /*@C 852 SNESLineSearchSet - Sets the line search routine to be used 853 by the method SNESLS. 854 855 Input Parameters: 856 + snes - nonlinear context obtained from SNESCreate() 857 . lsctx - optional user-defined context for use by line search 858 - func - pointer to int function 859 860 Collective on SNES 861 862 Available Routines: 863 + SNESLineSearchCubic() - default line search 864 . SNESLineSearchQuadratic() - quadratic line search 865 . SNESLineSearchNo() - the full Newton step (actually not a line search) 866 - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 867 868 Options Database Keys: 869 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 870 . -snes_ls_alpha <alpha> - Sets alpha 871 . -snes_ls_maxstep <maxstep> - Sets maximum step the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 872 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 873 874 Calling sequence of func: 875 .vb 876 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 877 .ve 878 879 Input parameters for func: 880 + snes - nonlinear context 881 . lsctx - optional user-defined context for line search 882 . x - current iterate 883 . f - residual evaluated at x 884 . y - search direction 885 . w - work vector 886 - fnorm - 2-norm of f 887 888 Output parameters for func: 889 + g - residual evaluated at new iterate y 890 . w - new iterate 891 . gnorm - 2-norm of g 892 . ynorm - 2-norm of search length 893 - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 894 895 Level: advanced 896 897 .keywords: SNES, nonlinear, set, line search, routine 898 899 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 900 SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 901 @*/ 902 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 903 { 904 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 905 906 PetscFunctionBegin; 907 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 908 if (f) { 909 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 910 } 911 PetscFunctionReturn(0); 912 } 913 914 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 915 /* -------------------------------------------------------------------------- */ 916 EXTERN_C_BEGIN 917 #undef __FUNCT__ 918 #define __FUNCT__ "SNESLineSearchSet_LS" 919 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 920 { 921 PetscFunctionBegin; 922 ((SNES_LS *)(snes->data))->LineSearch = func; 923 ((SNES_LS *)(snes->data))->lsP = lsctx; 924 PetscFunctionReturn(0); 925 } 926 EXTERN_C_END 927 /* -------------------------------------------------------------------------- */ 928 #undef __FUNCT__ 929 #define __FUNCT__ "SNESLineSearchSetPostCheck" 930 /*@C 931 SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 932 by the line search routine in the Newton-based method SNESLS. 933 934 Input Parameters: 935 + snes - nonlinear context obtained from SNESCreate() 936 . func - pointer to function 937 - checkctx - optional user-defined context for use by step checking routine 938 939 Collective on SNES 940 941 Calling sequence of func: 942 .vb 943 int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 944 .ve 945 where func returns an error code of 0 on success and a nonzero 946 on failure. 947 948 Input parameters for func: 949 + snes - nonlinear context 950 . checkctx - optional user-defined context for use by step checking routine 951 . x - previous iterate 952 . y - new search direction and length 953 - w - current candidate iterate 954 955 Output parameters for func: 956 + y - search direction (possibly changed) 957 . w - current iterate (possibly modified) 958 . changed_y - indicates search direction was changed by this routine 959 - changed_w - indicates current iterate was changed by this routine 960 961 Level: advanced 962 963 Notes: All line searches accept the new iterate computed by the line search checking routine. 964 965 Only one of changed_y and changed_w can be PETSC_TRUE 966 967 On input w = x + y 968 969 SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 970 to the checking routine, and then (3) compute the corresponding nonlinear 971 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 972 973 SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 974 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 975 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 976 were made to the candidate iterate in the checking routine (as indicated 977 by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 978 very costly, so use this feature with caution! 979 980 .keywords: SNES, nonlinear, set, line search check, step check, routine 981 982 .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 983 @*/ 984 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 985 { 986 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 987 988 PetscFunctionBegin; 989 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 990 if (f) { 991 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 992 } 993 PetscFunctionReturn(0); 994 } 995 #undef __FUNCT__ 996 #define __FUNCT__ "SNESLineSearchSetPreCheck" 997 /*@C 998 SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 999 before the line search is called. 1000 1001 Input Parameters: 1002 + snes - nonlinear context obtained from SNESCreate() 1003 . func - pointer to function 1004 - checkctx - optional user-defined context for use by step checking routine 1005 1006 Collective on SNES 1007 1008 Calling sequence of func: 1009 .vb 1010 int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 1011 .ve 1012 where func returns an error code of 0 on success and a nonzero 1013 on failure. 1014 1015 Input parameters for func: 1016 + snes - nonlinear context 1017 . checkctx - optional user-defined context for use by step checking routine 1018 . x - previous iterate 1019 - y - new search direction and length 1020 1021 Output parameters for func: 1022 + y - search direction (possibly changed) 1023 - changed_y - indicates search direction was changed by this routine 1024 1025 Level: advanced 1026 1027 Notes: All line searches accept the new iterate computed by the line search checking routine. 1028 1029 .keywords: SNES, nonlinear, set, line search check, step check, routine 1030 1031 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1032 @*/ 1033 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 1034 { 1035 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 1036 1037 PetscFunctionBegin; 1038 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1039 if (f) { 1040 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1041 } 1042 PetscFunctionReturn(0); 1043 } 1044 1045 /* -------------------------------------------------------------------------- */ 1046 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 1047 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1048 EXTERN_C_BEGIN 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1051 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1052 { 1053 PetscFunctionBegin; 1054 ((SNES_LS *)(snes->data))->postcheckstep = func; 1055 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1056 PetscFunctionReturn(0); 1057 } 1058 EXTERN_C_END 1059 1060 EXTERN_C_BEGIN 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 1063 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 1064 { 1065 PetscFunctionBegin; 1066 ((SNES_LS *)(snes->data))->precheckstep = func; 1067 ((SNES_LS *)(snes->data))->precheck = checkctx; 1068 PetscFunctionReturn(0); 1069 } 1070 EXTERN_C_END 1071 1072 /* 1073 SNESView_LS - Prints info from the SNESLS data structure. 1074 1075 Input Parameters: 1076 . SNES - the SNES context 1077 . viewer - visualization context 1078 1079 Application Interface Routine: SNESView() 1080 */ 1081 #undef __FUNCT__ 1082 #define __FUNCT__ "SNESView_LS" 1083 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1084 { 1085 SNES_LS *ls = (SNES_LS *)snes->data; 1086 const char *cstr; 1087 PetscErrorCode ierr; 1088 PetscTruth iascii; 1089 1090 PetscFunctionBegin; 1091 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1092 if (iascii) { 1093 if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 1094 else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 1095 else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 1096 else cstr = "unknown"; 1097 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1098 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr); 1099 } else { 1100 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1101 } 1102 PetscFunctionReturn(0); 1103 } 1104 /* -------------------------------------------------------------------------- */ 1105 /* 1106 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1107 1108 Input Parameter: 1109 . snes - the SNES context 1110 1111 Application Interface Routine: SNESSetFromOptions() 1112 */ 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "SNESSetFromOptions_LS" 1115 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1116 { 1117 SNES_LS *ls = (SNES_LS *)snes->data; 1118 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1119 PetscErrorCode ierr; 1120 PetscInt indx; 1121 PetscTruth flg; 1122 1123 PetscFunctionBegin; 1124 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1125 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1126 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1127 ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 1128 1129 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1130 if (flg) { 1131 switch (indx) { 1132 case 0: 1133 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1134 break; 1135 case 1: 1136 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1137 break; 1138 case 2: 1139 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1140 break; 1141 case 3: 1142 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1143 break; 1144 } 1145 } 1146 ierr = PetscOptionsTail();CHKERRQ(ierr); 1147 PetscFunctionReturn(0); 1148 } 1149 /* -------------------------------------------------------------------------- */ 1150 /*MC 1151 SNESLS - Newton based nonlinear solver that uses a line search 1152 1153 Options Database: 1154 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1155 . -snes_ls_alpha <alpha> - Sets alpha 1156 . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 1157 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1158 1159 Notes: This is the default nonlinear solver in SNES 1160 1161 Level: beginner 1162 1163 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1164 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1165 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1166 1167 M*/ 1168 EXTERN_C_BEGIN 1169 #undef __FUNCT__ 1170 #define __FUNCT__ "SNESCreate_LS" 1171 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 1172 { 1173 PetscErrorCode ierr; 1174 SNES_LS *neP; 1175 1176 PetscFunctionBegin; 1177 snes->ops->setup = SNESSetUp_LS; 1178 snes->ops->solve = SNESSolve_LS; 1179 snes->ops->destroy = SNESDestroy_LS; 1180 snes->ops->setfromoptions = SNESSetFromOptions_LS; 1181 snes->ops->view = SNESView_LS; 1182 1183 ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 1184 snes->data = (void*)neP; 1185 neP->alpha = 1.e-4; 1186 neP->maxstep = 1.e8; 1187 neP->minlambda = 1.e-12; 1188 neP->LineSearch = SNESLineSearchCubic; 1189 neP->lsP = PETSC_NULL; 1190 neP->postcheckstep = PETSC_NULL; 1191 neP->postcheck = PETSC_NULL; 1192 neP->precheckstep = PETSC_NULL; 1193 neP->precheck = PETSC_NULL; 1194 1195 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C", 1196 "SNESLineSearchSet_LS", 1197 SNESLineSearchSet_LS);CHKERRQ(ierr); 1198 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C", 1199 "SNESLineSearchSetPostCheck_LS", 1200 SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1201 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C", 1202 "SNESLineSearchSetPreCheck_LS", 1203 SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 1204 1205 PetscFunctionReturn(0); 1206 } 1207 EXTERN_C_END 1208 1209 1210 1211