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