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 != 0.0) { 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.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 . fnorm - 2-norm of f 344 - xnorm - norm of x if known, otherwise 0 345 346 Output Parameters: 347 + g - residual evaluated at new iterate y 348 . w - new iterate 349 . gnorm - 2-norm of g 350 . ynorm - 2-norm of search length 351 - flag - PETSC_TRUE on success, PETSC_FALSE on failure 352 353 Options Database Key: 354 . -snes_ls basic - Activates SNESLineSearchNo() 355 356 Level: advanced 357 358 .keywords: SNES, nonlinear, line search, cubic 359 360 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 361 SNESLineSearchSet(), SNESLineSearchNoNorms() 362 @*/ 363 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) 364 { 365 PetscErrorCode ierr; 366 SNES_LS *neP = (SNES_LS*)snes->data; 367 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 368 369 PetscFunctionBegin; 370 *flag = PETSC_TRUE; 371 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 372 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 373 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 374 if (neP->postcheckstep) { 375 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 376 } 377 if (changed_y) { 378 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 379 } 380 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 381 if (!snes->domainerror) { 382 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 383 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 384 } 385 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 386 PetscFunctionReturn(0); 387 } 388 /* -------------------------------------------------------------------------- */ 389 #undef __FUNCT__ 390 #define __FUNCT__ "SNESLineSearchNoNorms" 391 392 /*@C 393 SNESLineSearchNoNorms - This routine is not a line search at 394 all; it simply uses the full Newton step. This version does not 395 even compute the norm of the function or search direction; this 396 is intended only when you know the full step is fine and are 397 not checking for convergence of the nonlinear iteration (for 398 example, you are running always for a fixed number of Newton steps). 399 400 Collective on SNES and Vec 401 402 Input Parameters: 403 + snes - nonlinear context 404 . lsctx - optional context for line search (not used here) 405 . x - current iterate 406 . f - residual evaluated at x 407 . y - search direction 408 . w - work vector 409 . fnorm - 2-norm of f 410 - xnorm - norm of x if known, otherwise 0 411 412 Output Parameters: 413 + g - residual evaluated at new iterate y 414 . w - new iterate 415 . gnorm - not changed 416 . ynorm - not changed 417 - flag - set to PETSC_TRUE indicating a successful line search 418 419 Options Database Key: 420 . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 421 422 Notes: 423 SNESLineSearchNoNorms() must be used in conjunction with 424 either the options 425 $ -snes_no_convergence_test -snes_max_it <its> 426 or alternatively a user-defined custom test set via 427 SNESSetConvergenceTest(); or a -snes_max_it of 1, 428 otherwise, the SNES solver will generate an error. 429 430 During the final iteration this will not evaluate the function at 431 the solution point. This is to save a function evaluation while 432 using pseudo-timestepping. 433 434 The residual norms printed by monitoring routines such as 435 SNESMonitorDefault() (as activated via -snes_monitor) will not be 436 correct, since they are not computed. 437 438 Level: advanced 439 440 .keywords: SNES, nonlinear, line search, cubic 441 442 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 443 SNESLineSearchSet(), SNESLineSearchNo() 444 @*/ 445 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) 446 { 447 PetscErrorCode ierr; 448 SNES_LS *neP = (SNES_LS*)snes->data; 449 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 450 451 PetscFunctionBegin; 452 *flag = PETSC_TRUE; 453 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 454 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 455 if (neP->postcheckstep) { 456 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 457 } 458 if (changed_y) { 459 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 460 } 461 462 /* don't evaluate function the last time through */ 463 if (snes->iter < snes->max_its-1) { 464 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 465 } 466 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 467 PetscFunctionReturn(0); 468 } 469 /* -------------------------------------------------------------------------- */ 470 #undef __FUNCT__ 471 #define __FUNCT__ "SNESLineSearchCubic" 472 /*@C 473 SNESLineSearchCubic - Performs a cubic line search (default line search method). 474 475 Collective on SNES 476 477 Input Parameters: 478 + snes - nonlinear context 479 . lsctx - optional context for line search (not used here) 480 . x - current iterate 481 . f - residual evaluated at x 482 . y - search direction 483 . w - work vector 484 . fnorm - 2-norm of f 485 - xnorm - norm of x if known, otherwise 0 486 487 Output Parameters: 488 + g - residual evaluated at new iterate y 489 . w - new iterate 490 . gnorm - 2-norm of g 491 . ynorm - 2-norm of search length 492 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 493 494 Options Database Key: 495 + -snes_ls cubic - Activates SNESLineSearchCubic() 496 . -snes_ls_alpha <alpha> - Sets alpha 497 . -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) 498 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 499 500 501 Notes: 502 This line search is taken from "Numerical Methods for Unconstrained 503 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 504 505 Level: advanced 506 507 .keywords: SNES, nonlinear, line search, cubic 508 509 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 510 @*/ 511 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) 512 { 513 /* 514 Note that for line search purposes we work with with the related 515 minimization problem: 516 min z(x): R^n -> R, 517 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 518 */ 519 520 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 521 PetscReal minlambda,lambda,lambdatemp; 522 #if defined(PETSC_USE_COMPLEX) 523 PetscScalar cinitslope; 524 #endif 525 PetscErrorCode ierr; 526 PetscInt count; 527 SNES_LS *neP = (SNES_LS*)snes->data; 528 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 529 530 PetscFunctionBegin; 531 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 532 *flag = PETSC_TRUE; 533 534 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 535 if (*ynorm == 0.0) { 536 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 537 *gnorm = fnorm; 538 ierr = VecCopy(x,w);CHKERRQ(ierr); 539 ierr = VecCopy(f,g);CHKERRQ(ierr); 540 *flag = PETSC_FALSE; 541 goto theend1; 542 } 543 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 544 ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 545 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 546 *ynorm = neP->maxstep; 547 } 548 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 549 minlambda = neP->minlambda/rellength; 550 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 551 #if defined(PETSC_USE_COMPLEX) 552 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 553 initslope = PetscRealPart(cinitslope); 554 #else 555 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 556 #endif 557 if (initslope > 0.0) initslope = -initslope; 558 if (initslope == 0.0) initslope = -1.0; 559 560 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 561 if (snes->nfuncs >= snes->max_funcs) { 562 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 563 *flag = PETSC_FALSE; 564 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 565 goto theend1; 566 } 567 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 568 if (snes->domainerror) { 569 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 570 PetscFunctionReturn(0); 571 } 572 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 573 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 574 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 575 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 576 ierr = PetscInfo2(snes,"Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 577 goto theend1; 578 } 579 580 /* Fit points with quadratic */ 581 lambda = 1.0; 582 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 583 lambdaprev = lambda; 584 gnormprev = *gnorm; 585 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 586 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 587 else lambda = lambdatemp; 588 589 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 590 if (snes->nfuncs >= snes->max_funcs) { 591 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 592 *flag = PETSC_FALSE; 593 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 594 goto theend1; 595 } 596 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 597 if (snes->domainerror) { 598 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 599 PetscFunctionReturn(0); 600 } 601 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 602 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 603 ierr = PetscInfo1(snes,"gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 604 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 605 ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 606 goto theend1; 607 } 608 609 /* Fit points with cubic */ 610 count = 1; 611 while (PETSC_TRUE) { 612 if (lambda <= minlambda) { 613 ierr = PetscInfo1(snes,"Unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 614 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); 615 *flag = PETSC_FALSE; 616 break; 617 } 618 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 619 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 620 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 621 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 622 d = b*b - 3*a*initslope; 623 if (d < 0.0) d = 0.0; 624 if (a == 0.0) { 625 lambdatemp = -initslope/(2.0*b); 626 } else { 627 lambdatemp = (-b + sqrt(d))/(3.0*a); 628 } 629 lambdaprev = lambda; 630 gnormprev = *gnorm; 631 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 632 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 633 else lambda = lambdatemp; 634 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 635 if (snes->nfuncs >= snes->max_funcs) { 636 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 637 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); 638 *flag = PETSC_FALSE; 639 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 640 break; 641 } 642 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 643 if (snes->domainerror) { 644 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 645 PetscFunctionReturn(0); 646 } 647 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 648 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 649 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 650 ierr = PetscInfo2(snes,"Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 651 break; 652 } else { 653 ierr = PetscInfo2(snes,"Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 654 } 655 count++; 656 } 657 theend1: 658 /* Optional user-defined check for line search step validity */ 659 if (neP->postcheckstep && *flag) { 660 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 661 if (changed_y) { 662 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 663 } 664 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 665 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 666 if (snes->domainerror) { 667 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 668 PetscFunctionReturn(0); 669 } 670 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 671 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 672 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 673 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 674 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 675 } 676 } 677 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 678 PetscFunctionReturn(0); 679 } 680 /* -------------------------------------------------------------------------- */ 681 #undef __FUNCT__ 682 #define __FUNCT__ "SNESLineSearchQuadratic" 683 /*@C 684 SNESLineSearchQuadratic - Performs a quadratic line search. 685 686 Collective on SNES and Vec 687 688 Input Parameters: 689 + snes - the SNES context 690 . lsctx - optional context for line search (not used here) 691 . x - current iterate 692 . f - residual evaluated at x 693 . y - search direction 694 . w - work vector 695 . fnorm - 2-norm of f 696 - xnorm - norm of x if known, otherwise 0 697 698 Output Parameters: 699 + g - residual evaluated at new iterate w 700 . w - new iterate (x + lambda*y) 701 . gnorm - 2-norm of g 702 . ynorm - 2-norm of search length 703 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 704 705 Options Database Keys: 706 + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 707 . -snes_ls_alpha <alpha> - Sets alpha 708 . -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) 709 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 710 711 Notes: 712 Use SNESLineSearchSet() to set this routine within the SNESLS method. 713 714 Level: advanced 715 716 .keywords: SNES, nonlinear, quadratic, line search 717 718 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 719 @*/ 720 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) 721 { 722 /* 723 Note that for line search purposes we work with with the related 724 minimization problem: 725 min z(x): R^n -> R, 726 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 727 */ 728 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 729 #if defined(PETSC_USE_COMPLEX) 730 PetscScalar cinitslope; 731 #endif 732 PetscErrorCode ierr; 733 PetscInt count; 734 SNES_LS *neP = (SNES_LS*)snes->data; 735 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 736 737 PetscFunctionBegin; 738 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 739 *flag = PETSC_TRUE; 740 741 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 742 if (*ynorm == 0.0) { 743 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 744 *gnorm = fnorm; 745 ierr = VecCopy(x,w);CHKERRQ(ierr); 746 ierr = VecCopy(f,g);CHKERRQ(ierr); 747 *flag = PETSC_FALSE; 748 goto theend2; 749 } 750 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 751 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 752 *ynorm = neP->maxstep; 753 } 754 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 755 minlambda = neP->minlambda/rellength; 756 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 757 #if defined(PETSC_USE_COMPLEX) 758 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 759 initslope = PetscRealPart(cinitslope); 760 #else 761 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 762 #endif 763 if (initslope > 0.0) initslope = -initslope; 764 if (initslope == 0.0) initslope = -1.0; 765 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 766 767 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 768 if (snes->nfuncs >= snes->max_funcs) { 769 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 770 *flag = PETSC_FALSE; 771 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 772 goto theend2; 773 } 774 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 775 if (snes->domainerror) { 776 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 777 PetscFunctionReturn(0); 778 } 779 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 780 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 781 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 782 ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 783 goto theend2; 784 } 785 786 /* Fit points with quadratic */ 787 lambda = 1.0; 788 count = 1; 789 while (PETSC_TRUE) { 790 if (lambda <= minlambda) { /* bad luck; use full step */ 791 ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 792 ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 793 ierr = VecCopy(x,w);CHKERRQ(ierr); 794 *flag = PETSC_FALSE; 795 break; 796 } 797 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 798 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 799 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 800 else lambda = lambdatemp; 801 802 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 803 if (snes->nfuncs >= snes->max_funcs) { 804 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 805 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); 806 *flag = PETSC_FALSE; 807 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 808 break; 809 } 810 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 811 if (snes->domainerror) { 812 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 813 PetscFunctionReturn(0); 814 } 815 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 816 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 817 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 818 ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 819 break; 820 } 821 count++; 822 } 823 theend2: 824 /* Optional user-defined check for line search step validity */ 825 if (neP->postcheckstep) { 826 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 827 if (changed_y) { 828 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 829 } 830 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 831 ierr = SNESComputeFunction(snes,w,g); 832 if (snes->domainerror) { 833 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 834 PetscFunctionReturn(0); 835 } 836 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 837 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 838 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 839 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 840 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 841 } 842 } 843 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 844 PetscFunctionReturn(0); 845 } 846 847 /* -------------------------------------------------------------------------- */ 848 #undef __FUNCT__ 849 #define __FUNCT__ "SNESLineSearchSet" 850 /*@C 851 SNESLineSearchSet - Sets the line search routine to be used 852 by the method SNESLS. 853 854 Input Parameters: 855 + snes - nonlinear context obtained from SNESCreate() 856 . lsctx - optional user-defined context for use by line search 857 - func - pointer to int function 858 859 Collective on SNES 860 861 Available Routines: 862 + SNESLineSearchCubic() - default line search 863 . SNESLineSearchQuadratic() - quadratic line search 864 . SNESLineSearchNo() - the full Newton step (actually not a line search) 865 - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 866 867 Options Database Keys: 868 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 869 . -snes_ls_alpha <alpha> - Sets alpha 870 . -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) 871 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 872 873 Calling sequence of func: 874 .vb 875 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) 876 .ve 877 878 Input parameters for func: 879 + snes - nonlinear context 880 . lsctx - optional user-defined context for line search 881 . x - current iterate 882 . f - residual evaluated at x 883 . y - search direction 884 - fnorm - 2-norm of f 885 886 Output parameters for func: 887 + g - residual evaluated at new iterate y 888 . w - new iterate 889 . gnorm - 2-norm of g 890 . ynorm - 2-norm of search length 891 - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 892 893 Level: advanced 894 895 .keywords: SNES, nonlinear, set, line search, routine 896 897 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 898 SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 899 @*/ 900 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 901 { 902 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 903 904 PetscFunctionBegin; 905 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 906 if (f) { 907 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 908 } 909 PetscFunctionReturn(0); 910 } 911 912 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*/ 913 /* -------------------------------------------------------------------------- */ 914 EXTERN_C_BEGIN 915 #undef __FUNCT__ 916 #define __FUNCT__ "SNESLineSearchSet_LS" 917 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 918 { 919 PetscFunctionBegin; 920 ((SNES_LS *)(snes->data))->LineSearch = func; 921 ((SNES_LS *)(snes->data))->lsP = lsctx; 922 PetscFunctionReturn(0); 923 } 924 EXTERN_C_END 925 /* -------------------------------------------------------------------------- */ 926 #undef __FUNCT__ 927 #define __FUNCT__ "SNESLineSearchSetPostCheck" 928 /*@C 929 SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 930 by the line search routine in the Newton-based method SNESLS. 931 932 Input Parameters: 933 + snes - nonlinear context obtained from SNESCreate() 934 . func - pointer to function 935 - checkctx - optional user-defined context for use by step checking routine 936 937 Collective on SNES 938 939 Calling sequence of func: 940 .vb 941 int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 942 .ve 943 where func returns an error code of 0 on success and a nonzero 944 on failure. 945 946 Input parameters for func: 947 + snes - nonlinear context 948 . checkctx - optional user-defined context for use by step checking routine 949 . x - previous iterate 950 . y - new search direction and length 951 - w - current candidate iterate 952 953 Output parameters for func: 954 + y - search direction (possibly changed) 955 . w - current iterate (possibly modified) 956 . changed_y - indicates search direction was changed by this routine 957 - changed_w - indicates current iterate was changed by this routine 958 959 Level: advanced 960 961 Notes: All line searches accept the new iterate computed by the line search checking routine. 962 963 Only one of changed_y and changed_w can be PETSC_TRUE 964 965 On input w = x + y 966 967 SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 968 to the checking routine, and then (3) compute the corresponding nonlinear 969 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 970 971 SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 972 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 973 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 974 were made to the candidate iterate in the checking routine (as indicated 975 by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 976 very costly, so use this feature with caution! 977 978 .keywords: SNES, nonlinear, set, line search check, step check, routine 979 980 .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 981 @*/ 982 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 983 { 984 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 985 986 PetscFunctionBegin; 987 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 988 if (f) { 989 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 990 } 991 PetscFunctionReturn(0); 992 } 993 #undef __FUNCT__ 994 #define __FUNCT__ "SNESLineSearchSetPreCheck" 995 /*@C 996 SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 997 before the line search is called. 998 999 Input Parameters: 1000 + snes - nonlinear context obtained from SNESCreate() 1001 . func - pointer to function 1002 - checkctx - optional user-defined context for use by step checking routine 1003 1004 Collective on SNES 1005 1006 Calling sequence of func: 1007 .vb 1008 int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 1009 .ve 1010 where func returns an error code of 0 on success and a nonzero 1011 on failure. 1012 1013 Input parameters for func: 1014 + snes - nonlinear context 1015 . checkctx - optional user-defined context for use by step checking routine 1016 . x - previous iterate 1017 - y - new search direction and length 1018 1019 Output parameters for func: 1020 + y - search direction (possibly changed) 1021 - changed_y - indicates search direction was changed by this routine 1022 1023 Level: advanced 1024 1025 Notes: All line searches accept the new iterate computed by the line search checking routine. 1026 1027 .keywords: SNES, nonlinear, set, line search check, step check, routine 1028 1029 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1030 @*/ 1031 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 1032 { 1033 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 1034 1035 PetscFunctionBegin; 1036 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1037 if (f) { 1038 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1039 } 1040 PetscFunctionReturn(0); 1041 } 1042 1043 /* -------------------------------------------------------------------------- */ 1044 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 1045 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1046 EXTERN_C_BEGIN 1047 #undef __FUNCT__ 1048 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1049 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1050 { 1051 PetscFunctionBegin; 1052 ((SNES_LS *)(snes->data))->postcheckstep = func; 1053 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1054 PetscFunctionReturn(0); 1055 } 1056 EXTERN_C_END 1057 1058 EXTERN_C_BEGIN 1059 #undef __FUNCT__ 1060 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 1061 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 1062 { 1063 PetscFunctionBegin; 1064 ((SNES_LS *)(snes->data))->precheckstep = func; 1065 ((SNES_LS *)(snes->data))->precheck = checkctx; 1066 PetscFunctionReturn(0); 1067 } 1068 EXTERN_C_END 1069 1070 /* 1071 SNESView_LS - Prints info from the SNESLS data structure. 1072 1073 Input Parameters: 1074 . SNES - the SNES context 1075 . viewer - visualization context 1076 1077 Application Interface Routine: SNESView() 1078 */ 1079 #undef __FUNCT__ 1080 #define __FUNCT__ "SNESView_LS" 1081 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1082 { 1083 SNES_LS *ls = (SNES_LS *)snes->data; 1084 const char *cstr; 1085 PetscErrorCode ierr; 1086 PetscTruth iascii; 1087 1088 PetscFunctionBegin; 1089 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1090 if (iascii) { 1091 if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 1092 else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 1093 else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 1094 else cstr = "unknown"; 1095 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1096 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr); 1097 } else { 1098 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1099 } 1100 PetscFunctionReturn(0); 1101 } 1102 /* -------------------------------------------------------------------------- */ 1103 /* 1104 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1105 1106 Input Parameter: 1107 . snes - the SNES context 1108 1109 Application Interface Routine: SNESSetFromOptions() 1110 */ 1111 #undef __FUNCT__ 1112 #define __FUNCT__ "SNESSetFromOptions_LS" 1113 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1114 { 1115 SNES_LS *ls = (SNES_LS *)snes->data; 1116 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1117 PetscErrorCode ierr; 1118 PetscInt indx; 1119 PetscTruth flg; 1120 1121 PetscFunctionBegin; 1122 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1123 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1124 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1125 ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 1126 1127 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1128 if (flg) { 1129 switch (indx) { 1130 case 0: 1131 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1132 break; 1133 case 1: 1134 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1135 break; 1136 case 2: 1137 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1138 break; 1139 case 3: 1140 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1141 break; 1142 } 1143 } 1144 ierr = PetscOptionsTail();CHKERRQ(ierr); 1145 PetscFunctionReturn(0); 1146 } 1147 /* -------------------------------------------------------------------------- */ 1148 /*MC 1149 SNESLS - Newton based nonlinear solver that uses a line search 1150 1151 Options Database: 1152 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1153 . -snes_ls_alpha <alpha> - Sets alpha 1154 . -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) 1155 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1156 1157 Notes: This is the default nonlinear solver in SNES 1158 1159 Level: beginner 1160 1161 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1162 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1163 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1164 1165 M*/ 1166 EXTERN_C_BEGIN 1167 #undef __FUNCT__ 1168 #define __FUNCT__ "SNESCreate_LS" 1169 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 1170 { 1171 PetscErrorCode ierr; 1172 SNES_LS *neP; 1173 1174 PetscFunctionBegin; 1175 snes->ops->setup = SNESSetUp_LS; 1176 snes->ops->solve = SNESSolve_LS; 1177 snes->ops->destroy = SNESDestroy_LS; 1178 snes->ops->setfromoptions = SNESSetFromOptions_LS; 1179 snes->ops->view = SNESView_LS; 1180 1181 ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 1182 snes->data = (void*)neP; 1183 neP->alpha = 1.e-4; 1184 neP->maxstep = 1.e8; 1185 neP->minlambda = 1.e-12; 1186 neP->LineSearch = SNESLineSearchCubic; 1187 neP->lsP = PETSC_NULL; 1188 neP->postcheckstep = PETSC_NULL; 1189 neP->postcheck = PETSC_NULL; 1190 neP->precheckstep = PETSC_NULL; 1191 neP->precheck = PETSC_NULL; 1192 1193 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C", 1194 "SNESLineSearchSet_LS", 1195 SNESLineSearchSet_LS);CHKERRQ(ierr); 1196 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C", 1197 "SNESLineSearchSetPostCheck_LS", 1198 SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1199 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C", 1200 "SNESLineSearchSetPreCheck_LS", 1201 SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 1202 1203 PetscFunctionReturn(0); 1204 } 1205 EXTERN_C_END 1206 1207 1208 1209