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 static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin) 11 { 12 PetscReal a1; 13 PetscBool hastranspose; 14 Vec W; 15 PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 16 17 PetscFunctionBegin; 18 *ismin = PETSC_FALSE; 19 PetscCall(SNESGetObjective(snes, &objective, NULL)); 20 if (!objective) { 21 PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 22 PetscCall(VecDuplicate(F, &W)); 23 if (hastranspose) { 24 /* Compute || J^T F|| */ 25 PetscCall(MatMultTranspose(A, F, W)); 26 PetscCall(VecNorm(W, NORM_2, &a1)); 27 PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm))); 28 if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 29 } else { 30 Vec work; 31 PetscScalar result; 32 PetscReal wnorm; 33 34 PetscCall(VecSetRandom(W, NULL)); 35 PetscCall(VecNorm(W, NORM_2, &wnorm)); 36 PetscCall(VecDuplicate(W, &work)); 37 PetscCall(MatMult(A, W, work)); 38 PetscCall(VecDot(F, work, &result)); 39 PetscCall(VecDestroy(&work)); 40 a1 = PetscAbsScalar(result) / (fnorm * wnorm); 41 PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1)); 42 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 43 } 44 PetscCall(VecDestroy(&W)); 45 } 46 PetscFunctionReturn(PETSC_SUCCESS); 47 } 48 49 /* 50 Checks if J^T(F - J*X) = 0 51 */ 52 static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X) 53 { 54 PetscReal a1, a2; 55 PetscBool hastranspose; 56 PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 57 58 PetscFunctionBegin; 59 PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 60 PetscCall(SNESGetObjective(snes, &objective, NULL)); 61 if (hastranspose && !objective) { 62 Vec W1, W2; 63 64 PetscCall(VecDuplicate(F, &W1)); 65 PetscCall(VecDuplicate(F, &W2)); 66 PetscCall(MatMult(A, X, W1)); 67 PetscCall(VecAXPY(W1, -1.0, F)); 68 69 /* Compute || J^T W|| */ 70 PetscCall(MatMultTranspose(A, W1, W2)); 71 PetscCall(VecNorm(W1, NORM_2, &a1)); 72 PetscCall(VecNorm(W2, NORM_2, &a2)); 73 if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n", (double)(a2 / a1))); 74 PetscCall(VecDestroy(&W1)); 75 PetscCall(VecDestroy(&W2)); 76 } 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79 80 /* 81 This file implements a truncated Newton method with a line search, 82 for solving a system of nonlinear equations, using the KSP, Vec, 83 and Mat interfaces for linear solvers, vectors, and matrices, 84 respectively. 85 86 The following basic routines are required for each nonlinear solver: 87 SNESCreate_XXX() - Creates a nonlinear solver context 88 SNESSetFromOptions_XXX() - Sets runtime options 89 SNESSolve_XXX() - Solves the nonlinear system 90 SNESDestroy_XXX() - Destroys the nonlinear solver context 91 The suffix "_XXX" denotes a particular implementation, in this case 92 we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 93 systems of nonlinear equations with a line search (LS) method. 94 These routines are actually called via the common user interface 95 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 96 SNESDestroy(), so the application code interface remains identical 97 for all nonlinear solvers. 98 99 Another key routine is: 100 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 101 by setting data structures and options. The interface routine SNESSetUp() 102 is not usually called directly by the user, but instead is called by 103 SNESSolve() if necessary. 104 105 Additional basic routines are: 106 SNESView_XXX() - Prints details of runtime options that 107 have actually been used. 108 These are called by application codes via the interface routines 109 SNESView(). 110 111 The various types of solvers (preconditioners, Krylov subspace methods, 112 nonlinear solvers, timesteppers) are all organized similarly, so the 113 above description applies to these categories also. 114 115 */ 116 /* 117 SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 118 method with a line search. 119 120 Input Parameter: 121 . snes - the SNES context 122 123 Output Parameter: 124 . outits - number of iterations until termination 125 126 Application Interface Routine: SNESSolve() 127 128 */ 129 PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 130 { 131 PetscInt maxits, i, lits; 132 SNESLineSearchReason lssucceed; 133 PetscReal fnorm, xnorm, ynorm; 134 Vec Y, X, F; 135 SNESLineSearch linesearch; 136 SNESConvergedReason reason; 137 #if defined(PETSC_USE_INFO) 138 PetscReal gnorm; 139 #endif 140 141 PetscFunctionBegin; 142 PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 143 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->vec_sol_update; /* newton step */ 152 153 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 154 snes->iter = 0; 155 snes->norm = 0.0; 156 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 157 PetscCall(SNESGetLineSearch(snes, &linesearch)); 158 159 /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 160 if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 161 PetscCall(SNESApplyNPC(snes, X, NULL, F)); 162 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 163 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 164 snes->reason = SNES_DIVERGED_INNER; 165 PetscFunctionReturn(PETSC_SUCCESS); 166 } 167 168 PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 169 PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 170 } else { 171 if (!snes->vec_func_init_set) { 172 PetscCall(SNESComputeFunction(snes, X, F)); 173 } else snes->vec_func_init_set = PETSC_FALSE; 174 } 175 176 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 177 SNESCheckFunctionNorm(snes, fnorm); 178 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 179 snes->norm = fnorm; 180 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 181 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 182 183 /* test convergence */ 184 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 185 PetscCall(SNESMonitor(snes, 0, fnorm)); 186 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 187 188 for (i = 0; i < maxits; i++) { 189 /* Call general purpose update function */ 190 PetscTryTypeMethod(snes, update, snes->iter); 191 192 /* apply the nonlinear preconditioner */ 193 if (snes->npc) { 194 if (snes->npcside == PC_RIGHT) { 195 PetscCall(SNESSetInitialFunction(snes->npc, F)); 196 PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 197 PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 198 PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 199 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 200 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 201 snes->reason = SNES_DIVERGED_INNER; 202 PetscFunctionReturn(PETSC_SUCCESS); 203 } 204 PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 205 } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 206 PetscCall(SNESApplyNPC(snes, X, F, F)); 207 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 208 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 209 snes->reason = SNES_DIVERGED_INNER; 210 PetscFunctionReturn(PETSC_SUCCESS); 211 } 212 } 213 } 214 215 /* Solve J Y = F, where J is Jacobian matrix */ 216 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 217 SNESCheckJacobianDomainerror(snes); 218 PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 219 PetscCall(KSPSolve(snes->ksp, F, Y)); 220 SNESCheckKSPSolve(snes); 221 PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 222 PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 223 224 if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y)); 225 226 #if defined(PETSC_USE_INFO) 227 gnorm = fnorm; 228 #endif 229 /* Compute a (scaled) negative update in the line search routine: 230 X <- X - lambda*Y 231 and evaluate F = function(X) (depends on the line search). 232 */ 233 PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 234 PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed)); 235 PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 236 PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed)); 237 if (snes->reason) break; 238 SNESCheckFunctionNorm(snes, fnorm); 239 if (lssucceed) { 240 if (snes->stol * xnorm > ynorm) { 241 snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 242 PetscFunctionReturn(PETSC_SUCCESS); 243 } 244 if (++snes->numFailures >= snes->maxFailures) { 245 PetscBool ismin; 246 snes->reason = SNES_DIVERGED_LINE_SEARCH; 247 PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin)); 248 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 249 if (snes->errorifnotconverged && snes->reason) { 250 PetscViewer monitor; 251 PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 252 PetscCheck(monitor, PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due to %s. Suggest running with -snes_linesearch_monitor", SNESConvergedReasons[snes->reason]); 253 SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]); 254 } 255 break; 256 } 257 } 258 /* Monitor convergence */ 259 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 260 snes->iter = i + 1; 261 snes->norm = fnorm; 262 snes->ynorm = ynorm; 263 snes->xnorm = xnorm; 264 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 265 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 266 /* Test for convergence */ 267 PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 268 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 269 if (snes->reason) break; 270 } 271 PetscFunctionReturn(PETSC_SUCCESS); 272 } 273 274 /* 275 SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 276 of the SNESNEWTONLS nonlinear solver. 277 278 Input Parameter: 279 . snes - the SNES context 280 . x - the solution vector 281 282 Application Interface Routine: SNESSetUp() 283 284 */ 285 PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 286 { 287 PetscFunctionBegin; 288 PetscCall(SNESSetUpMatrices(snes)); 289 if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 290 PetscFunctionReturn(PETSC_SUCCESS); 291 } 292 293 PetscErrorCode SNESReset_NEWTONLS(SNES snes) 294 { 295 PetscFunctionBegin; 296 PetscFunctionReturn(PETSC_SUCCESS); 297 } 298 299 /* 300 SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 301 with SNESCreate_NEWTONLS(). 302 303 Input Parameter: 304 . snes - the SNES context 305 306 Application Interface Routine: SNESDestroy() 307 */ 308 PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 309 { 310 PetscFunctionBegin; 311 PetscCall(SNESReset_NEWTONLS(snes)); 312 PetscCall(PetscFree(snes->data)); 313 PetscFunctionReturn(PETSC_SUCCESS); 314 } 315 316 /* 317 SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 318 319 Input Parameters: 320 . SNES - the SNES context 321 . viewer - visualization context 322 323 Application Interface Routine: SNESView() 324 */ 325 static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) 326 { 327 PetscBool iascii; 328 329 PetscFunctionBegin; 330 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 331 if (iascii) { } 332 PetscFunctionReturn(PETSC_SUCCESS); 333 } 334 335 /* 336 SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 337 338 Input Parameter: 339 . snes - the SNES context 340 341 Application Interface Routine: SNESSetFromOptions() 342 */ 343 static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) 344 { 345 PetscFunctionBegin; 346 PetscFunctionReturn(PETSC_SUCCESS); 347 } 348 349 /*MC 350 SNESNEWTONLS - Newton based nonlinear solver that uses a line search 351 352 Options Database Keys: 353 + -snes_linesearch_type <bt> - bt,basic. Select line search type 354 . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 355 . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`) 356 . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 357 . -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) 358 . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 359 . -snes_linesearch_monitor - print information about progress of line searches 360 - -snes_linesearch_damping - damping factor used for basic line search 361 362 Note: 363 This is the default nonlinear solver in `SNES` 364 365 Level: beginner 366 367 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()` 368 `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()` 369 M*/ 370 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 371 { 372 SNES_NEWTONLS *neP; 373 SNESLineSearch linesearch; 374 375 PetscFunctionBegin; 376 snes->ops->setup = SNESSetUp_NEWTONLS; 377 snes->ops->solve = SNESSolve_NEWTONLS; 378 snes->ops->destroy = SNESDestroy_NEWTONLS; 379 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 380 snes->ops->view = SNESView_NEWTONLS; 381 snes->ops->reset = SNESReset_NEWTONLS; 382 383 snes->npcside = PC_RIGHT; 384 snes->usesksp = PETSC_TRUE; 385 snes->usesnpc = PETSC_TRUE; 386 387 PetscCall(SNESGetLineSearch(snes, &linesearch)); 388 if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 389 390 snes->alwayscomputesfinalresidual = PETSC_TRUE; 391 392 PetscCall(PetscNew(&neP)); 393 snes->data = (void *)neP; 394 PetscFunctionReturn(PETSC_SUCCESS); 395 } 396