1 2 #include <../src/snes/impls/ls/lsimpl.h> 3 4 /* 5 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 6 || 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 7 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 8 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 9 */ 10 #undef __FUNCT__ 11 #define __FUNCT__ "SNESNEWTONLSCheckLocalMin_Private" 12 PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 13 { 14 PetscReal a1; 15 PetscErrorCode ierr; 16 PetscBool hastranspose; 17 18 PetscFunctionBegin; 19 *ismin = PETSC_FALSE; 20 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 21 if (hastranspose) { 22 /* Compute || J^T F|| */ 23 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 24 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 25 ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 26 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 27 } else { 28 Vec work; 29 PetscScalar result; 30 PetscReal wnorm; 31 32 ierr = VecSetRandom(W,NULL);CHKERRQ(ierr); 33 ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 34 ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 35 ierr = MatMult(A,W,work);CHKERRQ(ierr); 36 ierr = VecDot(F,work,&result);CHKERRQ(ierr); 37 ierr = VecDestroy(&work);CHKERRQ(ierr); 38 a1 = PetscAbsScalar(result)/(fnorm*wnorm); 39 ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 40 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 41 } 42 PetscFunctionReturn(0); 43 } 44 45 /* 46 Checks if J^T(F - J*X) = 0 47 */ 48 #undef __FUNCT__ 49 #define __FUNCT__ "SNESNEWTONLSCheckResidual_Private" 50 PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 51 { 52 PetscReal a1,a2; 53 PetscErrorCode ierr; 54 PetscBool hastranspose; 55 56 PetscFunctionBegin; 57 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 58 if (hastranspose) { 59 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 60 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 61 62 /* Compute || J^T W|| */ 63 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 64 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 65 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 66 if (a1 != 0.0) { 67 ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 68 } 69 } 70 PetscFunctionReturn(0); 71 } 72 73 /* -------------------------------------------------------------------- 74 75 This file implements a truncated Newton method with a line search, 76 for solving a system of nonlinear equations, using the KSP, Vec, 77 and Mat interfaces for linear solvers, vectors, and matrices, 78 respectively. 79 80 The following basic routines are required for each nonlinear solver: 81 SNESCreate_XXX() - Creates a nonlinear solver context 82 SNESSetFromOptions_XXX() - Sets runtime options 83 SNESSolve_XXX() - Solves the nonlinear system 84 SNESDestroy_XXX() - Destroys the nonlinear solver context 85 The suffix "_XXX" denotes a particular implementation, in this case 86 we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 87 systems of nonlinear equations with a line search (LS) method. 88 These routines are actually called via the common user interface 89 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 90 SNESDestroy(), so the application code interface remains identical 91 for all nonlinear solvers. 92 93 Another key routine is: 94 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 95 by setting data structures and options. The interface routine SNESSetUp() 96 is not usually called directly by the user, but instead is called by 97 SNESSolve() if necessary. 98 99 Additional basic routines are: 100 SNESView_XXX() - Prints details of runtime options that 101 have actually been used. 102 These are called by application codes via the interface routines 103 SNESView(). 104 105 The various types of solvers (preconditioners, Krylov subspace methods, 106 nonlinear solvers, timesteppers) are all organized similarly, so the 107 above description applies to these categories also. 108 109 -------------------------------------------------------------------- */ 110 /* 111 SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 112 method with a line search. 113 114 Input Parameters: 115 . snes - the SNES context 116 117 Output Parameter: 118 . outits - number of iterations until termination 119 120 Application Interface Routine: SNESSolve() 121 122 Notes: 123 This implements essentially a truncated Newton method with a 124 line search. By default a cubic backtracking line search 125 is employed, as described in the text "Numerical Methods for 126 Unconstrained Optimization and Nonlinear Equations" by Dennis 127 and Schnabel. 128 */ 129 #undef __FUNCT__ 130 #define __FUNCT__ "SNESSolve_NEWTONLS" 131 PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 132 { 133 PetscErrorCode ierr; 134 PetscInt maxits,i,lits; 135 PetscBool lssucceed; 136 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 137 PetscReal fnorm,gnorm,xnorm,ynorm; 138 Vec Y,X,F,G,W,FPC; 139 KSPConvergedReason kspreason; 140 PetscBool domainerror; 141 SNESLineSearch linesearch; 142 SNESConvergedReason reason; 143 144 PetscFunctionBegin; 145 snes->numFailures = 0; 146 snes->numLinearSolveFailures = 0; 147 snes->reason = SNES_CONVERGED_ITERATING; 148 149 maxits = snes->max_its; /* maximum number of iterations */ 150 X = snes->vec_sol; /* solution vector */ 151 F = snes->vec_func; /* residual vector */ 152 Y = snes->vec_sol_update; /* newton step */ 153 G = snes->work[0]; 154 W = snes->work[1]; 155 156 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 157 snes->iter = 0; 158 snes->norm = 0.0; 159 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 160 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 161 162 /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 163 if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 164 ierr = SNESApplyPC(snes,X,PETSC_NULL,PETSC_NULL,F);CHKERRQ(ierr); 165 ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); 166 ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 167 } else { 168 if (!snes->vec_func_init_set) { 169 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 170 ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 171 if (domainerror) { 172 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 173 PetscFunctionReturn(0); 174 } 175 } else snes->vec_func_init_set = PETSC_FALSE; 176 177 if (!snes->norm_init_set) { 178 ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 179 ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 180 if (PetscIsInfOrNanReal(fnorm)) { 181 snes->reason = SNES_DIVERGED_FNORM_NAN; 182 PetscFunctionReturn(0); 183 } 184 } else { 185 fnorm = snes->norm_init; 186 snes->norm_init_set = PETSC_FALSE; 187 } 188 } 189 190 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 191 snes->norm = fnorm; 192 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 193 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 194 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 195 196 /* set parameter for default relative tolerance convergence test */ 197 snes->ttol = fnorm*snes->rtol; 198 /* test convergence */ 199 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 200 if (snes->reason) PetscFunctionReturn(0); 201 202 for (i=0; i<maxits; i++) { 203 204 /* Call general purpose update function */ 205 if (snes->ops->update) { 206 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 207 } 208 209 /* apply the nonlinear preconditioner */ 210 if (snes->pc) { 211 if (snes->pcside == PC_RIGHT) { 212 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 213 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 214 ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr); 215 ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr); 216 ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr); 217 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 218 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 219 snes->reason = SNES_DIVERGED_INNER; 220 PetscFunctionReturn(0); 221 } 222 ierr = SNESGetFunction(snes->pc, &FPC, NULL, NULL);CHKERRQ(ierr); 223 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 224 ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 225 } else if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 226 ierr = SNESApplyPC(snes,X,F,&fnorm,F);CHKERRQ(ierr); 227 } 228 } 229 230 /* Solve J Y = F, where J is Jacobian matrix */ 231 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 232 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 233 ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 234 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 235 if (kspreason < 0) { 236 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 237 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 238 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 239 break; 240 } 241 } 242 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 243 snes->linear_its += lits; 244 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 245 246 if (PetscLogPrintInfo) { 247 ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 248 } 249 250 /* Compute a (scaled) negative update in the line search routine: 251 X <- X - lambda*Y 252 and evaluate F = function(X) (depends on the line search). 253 */ 254 gnorm = fnorm; 255 ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 256 ierr = SNESLineSearchGetSuccess(linesearch, &lssucceed);CHKERRQ(ierr); 257 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 258 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)gnorm,(double)fnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 259 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 260 ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 261 if (domainerror) { 262 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 263 PetscFunctionReturn(0); 264 } 265 if (!lssucceed) { 266 if (snes->stol*xnorm > ynorm) { 267 snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 268 PetscFunctionReturn(0); 269 } 270 if (++snes->numFailures >= snes->maxFailures) { 271 PetscBool ismin; 272 snes->reason = SNES_DIVERGED_LINE_SEARCH; 273 ierr = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 274 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 275 break; 276 } 277 } 278 /* Monitor convergence */ 279 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 280 snes->iter = i+1; 281 snes->norm = fnorm; 282 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 283 ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 284 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 285 /* Test for convergence */ 286 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 287 if (snes->reason) break; 288 } 289 if (i == maxits) { 290 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 291 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 292 } 293 PetscFunctionReturn(0); 294 } 295 /* -------------------------------------------------------------------------- */ 296 /* 297 SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 298 of the SNESNEWTONLS nonlinear solver. 299 300 Input Parameter: 301 . snes - the SNES context 302 . x - the solution vector 303 304 Application Interface Routine: SNESSetUp() 305 306 Notes: 307 For basic use of the SNES solvers, the user need not explicitly call 308 SNESSetUp(), since these actions will automatically occur during 309 the call to SNESSolve(). 310 */ 311 #undef __FUNCT__ 312 #define __FUNCT__ "SNESSetUp_NEWTONLS" 313 PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 314 { 315 PetscErrorCode ierr; 316 317 PetscFunctionBegin; 318 ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr); 319 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 320 PetscFunctionReturn(0); 321 } 322 /* -------------------------------------------------------------------------- */ 323 324 #undef __FUNCT__ 325 #define __FUNCT__ "SNESReset_NEWTONLS" 326 PetscErrorCode SNESReset_NEWTONLS(SNES snes) 327 { 328 PetscFunctionBegin; 329 PetscFunctionReturn(0); 330 } 331 332 /* 333 SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 334 with SNESCreate_NEWTONLS(). 335 336 Input Parameter: 337 . snes - the SNES context 338 339 Application Interface Routine: SNESDestroy() 340 */ 341 #undef __FUNCT__ 342 #define __FUNCT__ "SNESDestroy_NEWTONLS" 343 PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 344 { 345 PetscErrorCode ierr; 346 347 PetscFunctionBegin; 348 ierr = SNESReset_NEWTONLS(snes);CHKERRQ(ierr); 349 ierr = PetscFree(snes->data);CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352 /* -------------------------------------------------------------------------- */ 353 354 /* 355 SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 356 357 Input Parameters: 358 . SNES - the SNES context 359 . viewer - visualization context 360 361 Application Interface Routine: SNESView() 362 */ 363 #undef __FUNCT__ 364 #define __FUNCT__ "SNESView_NEWTONLS" 365 static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer) 366 { 367 PetscErrorCode ierr; 368 PetscBool iascii; 369 370 PetscFunctionBegin; 371 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 372 if (iascii) { 373 } 374 PetscFunctionReturn(0); 375 } 376 377 /* -------------------------------------------------------------------------- */ 378 /* 379 SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 380 381 Input Parameter: 382 . snes - the SNES context 383 384 Application Interface Routine: SNESSetFromOptions() 385 */ 386 #undef __FUNCT__ 387 #define __FUNCT__ "SNESSetFromOptions_NEWTONLS" 388 static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes) 389 { 390 PetscErrorCode ierr; 391 SNESLineSearch linesearch; 392 393 PetscFunctionBegin; 394 ierr = PetscOptionsHead("SNESNEWTONLS options");CHKERRQ(ierr); 395 ierr = PetscOptionsTail();CHKERRQ(ierr); 396 /* set the default line search type */ 397 if (!snes->linesearch) { 398 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 399 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 400 } 401 PetscFunctionReturn(0); 402 } 403 404 /* -------------------------------------------------------------------------- */ 405 /*MC 406 SNESNEWTONLS - Newton based nonlinear solver that uses a line search 407 408 Options Database: 409 + -snes_linesearch_type <bt> - bt,basic. Select line search type 410 . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 411 . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch 412 . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 413 . -snes_linesearch_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) 414 . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 415 . -snes_linesearch_monitor - print information about progress of line searches 416 - -snes_linesearch_damping - damping factor used for basic line search 417 418 Notes: This is the default nonlinear solver in SNES 419 420 Level: beginner 421 422 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder() 423 SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms() 424 425 M*/ 426 #undef __FUNCT__ 427 #define __FUNCT__ "SNESCreate_NEWTONLS" 428 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 429 { 430 PetscErrorCode ierr; 431 SNES_NEWTONLS *neP; 432 433 PetscFunctionBegin; 434 snes->ops->setup = SNESSetUp_NEWTONLS; 435 snes->ops->solve = SNESSolve_NEWTONLS; 436 snes->ops->destroy = SNESDestroy_NEWTONLS; 437 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 438 snes->ops->view = SNESView_NEWTONLS; 439 snes->ops->reset = SNESReset_NEWTONLS; 440 441 snes->usesksp = PETSC_TRUE; 442 snes->usespc = PETSC_TRUE; 443 ierr = PetscNewLog(snes,SNES_NEWTONLS,&neP);CHKERRQ(ierr); 444 snes->data = (void*)neP; 445 PetscFunctionReturn(0); 446 } 447