1 #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/ 2 const char *const SNESNCGTypes[] = {"FR", "PRP", "HS", "DY", "CD", "SNESNCGType", "SNES_NCG_", NULL}; 3 4 static PetscErrorCode SNESReset_NCG(SNES snes) 5 { 6 PetscFunctionBegin; 7 PetscFunctionReturn(PETSC_SUCCESS); 8 } 9 10 static PetscErrorCode SNESDestroy_NCG(SNES snes) 11 { 12 PetscFunctionBegin; 13 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL)); 14 PetscCall(PetscFree(snes->data)); 15 PetscFunctionReturn(PETSC_SUCCESS); 16 } 17 18 static PetscErrorCode SNESSetUp_NCG(SNES snes) 19 { 20 PetscFunctionBegin; 21 PetscCall(SNESSetWorkVecs(snes, 2)); 22 PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 23 if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 24 PetscFunctionReturn(PETSC_SUCCESS); 25 } 26 27 static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 28 { 29 PetscScalar alpha, ptAp; 30 Vec X, Y, F, W; 31 SNES snes; 32 PetscReal *fnorm, *xnorm, *ynorm; 33 34 PetscFunctionBegin; 35 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 36 X = linesearch->vec_sol; 37 W = linesearch->vec_sol_new; 38 F = linesearch->vec_func; 39 Y = linesearch->vec_update; 40 fnorm = &linesearch->fnorm; 41 xnorm = &linesearch->xnorm; 42 ynorm = &linesearch->ynorm; 43 44 if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes)); 45 46 /* 47 The exact step size for unpreconditioned linear CG is just: 48 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 49 */ 50 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 51 PetscCall(VecDot(F, F, &alpha)); 52 PetscCall(MatMult(snes->jacobian, Y, W)); 53 PetscCall(VecDot(Y, W, &ptAp)); 54 alpha = alpha / ptAp; 55 PetscCall(VecAXPY(X, -alpha, Y)); 56 PetscCall(SNESComputeFunction(snes, X, F)); 57 58 PetscCall(VecNorm(F, NORM_2, fnorm)); 59 PetscCall(VecNorm(X, NORM_2, xnorm)); 60 PetscCall(VecNorm(Y, NORM_2, ynorm)); 61 PetscFunctionReturn(PETSC_SUCCESS); 62 } 63 64 /*MC 65 SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG` 66 67 This line search uses the length "as if" the problem is linear (that is what is computed by the linear CG method) using the Jacobian of the function. 68 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) where r (f) is the current residual (function value), p (y) is the current search direction. 69 70 Notes: 71 This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian 72 73 This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone. 74 75 Level: advanced 76 77 .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()` 78 M*/ 79 80 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 81 { 82 PetscFunctionBegin; 83 linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 84 linesearch->ops->destroy = NULL; 85 linesearch->ops->setfromoptions = NULL; 86 linesearch->ops->reset = NULL; 87 linesearch->ops->view = NULL; 88 linesearch->ops->setup = NULL; 89 PetscFunctionReturn(PETSC_SUCCESS); 90 } 91 92 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject) 93 { 94 SNES_NCG *ncg = (SNES_NCG *)snes->data; 95 PetscBool debug = PETSC_FALSE; 96 SNESNCGType ncgtype = ncg->type; 97 SNESLineSearch linesearch; 98 99 PetscFunctionBegin; 100 PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options"); 101 PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL)); 102 if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 103 PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL)); 104 PetscCall(SNESNCGSetType(snes, ncgtype)); 105 PetscOptionsHeadEnd(); 106 if (!snes->linesearch) { 107 PetscCall(SNESGetLineSearch(snes, &linesearch)); 108 if (!((PetscObject)linesearch)->type_name) { 109 if (!snes->npc) { 110 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP)); 111 } else { 112 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 113 } 114 } 115 } 116 PetscFunctionReturn(PETSC_SUCCESS); 117 } 118 119 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 120 { 121 SNES_NCG *ncg = (SNES_NCG *)snes->data; 122 PetscBool iascii; 123 124 PetscFunctionBegin; 125 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 126 if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type])); 127 PetscFunctionReturn(PETSC_SUCCESS); 128 } 129 130 /* 131 132 Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 133 134 This routine is not currently used. I don't know what its intended purpose is. 135 136 Note it has a hardwired differencing parameter of 1e-5 137 138 */ 139 PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar *ytJtf) 140 { 141 PetscScalar ftf, ftg, fty, h; 142 143 PetscFunctionBegin; 144 PetscCall(VecDot(F, F, &ftf)); 145 PetscCall(VecDot(F, Y, &fty)); 146 h = 1e-5 * fty / fty; 147 PetscCall(VecCopy(X, W)); 148 PetscCall(VecAXPY(W, -h, Y)); /* this is arbitrary */ 149 PetscCall(SNESComputeFunction(snes, W, G)); 150 PetscCall(VecDot(G, F, &ftg)); 151 *ytJtf = (ftg - ftf) / h; 152 PetscFunctionReturn(PETSC_SUCCESS); 153 } 154 155 /*@ 156 SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`. 157 158 Logically Collective 159 160 Input Parameters: 161 + snes - the iterative context 162 - btype - update type, see `SNESNCGType` 163 164 Options Database Key: 165 . -snes_ncg_type <prp,fr,hs,dy,cd> - strategy for selecting algorithm for computing beta 166 167 Level: intermediate 168 169 Notes: 170 `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions. 171 172 It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner, 173 that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`? 174 175 Developer Notes: 176 There should be a `SNESNCGSetType()` 177 178 .seealso: `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD` 179 @*/ 180 PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 181 { 182 PetscFunctionBegin; 183 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 184 PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype)); 185 PetscFunctionReturn(PETSC_SUCCESS); 186 } 187 188 static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 189 { 190 SNES_NCG *ncg = (SNES_NCG *)snes->data; 191 192 PetscFunctionBegin; 193 ncg->type = btype; 194 PetscFunctionReturn(PETSC_SUCCESS); 195 } 196 197 /* 198 SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 199 200 Input Parameter: 201 . snes - the `SNES` context 202 203 Output Parameter: 204 . outits - number of iterations until termination 205 206 Application Interface Routine: SNESSolve() 207 */ 208 static PetscErrorCode SNESSolve_NCG(SNES snes) 209 { 210 SNES_NCG *ncg = (SNES_NCG *)snes->data; 211 Vec X, dX, lX, F, dXold; 212 PetscReal fnorm, ynorm, xnorm, beta = 0.0; 213 PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 214 PetscInt maxits, i; 215 SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 216 SNESLineSearch linesearch; 217 SNESConvergedReason reason; 218 219 PetscFunctionBegin; 220 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); 221 222 PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 223 snes->reason = SNES_CONVERGED_ITERATING; 224 225 maxits = snes->max_its; /* maximum number of iterations */ 226 X = snes->vec_sol; /* X^n */ 227 dXold = snes->work[0]; /* The previous iterate of X */ 228 dX = snes->work[1]; /* the preconditioned direction */ 229 lX = snes->vec_sol_update; /* the conjugate direction */ 230 F = snes->vec_func; /* residual vector */ 231 232 PetscCall(SNESGetLineSearch(snes, &linesearch)); 233 234 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 235 snes->iter = 0; 236 snes->norm = 0.; 237 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 238 239 /* compute the initial function and preconditioned update dX */ 240 241 if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 242 PetscCall(SNESApplyNPC(snes, X, NULL, dX)); 243 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 244 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 245 snes->reason = SNES_DIVERGED_INNER; 246 PetscFunctionReturn(PETSC_SUCCESS); 247 } 248 PetscCall(VecCopy(dX, F)); 249 PetscCall(VecNorm(F, NORM_2, &fnorm)); 250 } else { 251 if (!snes->vec_func_init_set) { 252 PetscCall(SNESComputeFunction(snes, X, F)); 253 } else snes->vec_func_init_set = PETSC_FALSE; 254 255 /* convergence test */ 256 PetscCall(VecNorm(F, NORM_2, &fnorm)); 257 SNESCheckFunctionNorm(snes, fnorm); 258 PetscCall(VecCopy(F, dX)); 259 } 260 if (snes->npc) { 261 if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 262 PetscCall(SNESApplyNPC(snes, X, F, dX)); 263 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 264 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 265 snes->reason = SNES_DIVERGED_INNER; 266 PetscFunctionReturn(PETSC_SUCCESS); 267 } 268 } 269 } 270 PetscCall(VecCopy(dX, lX)); 271 PetscCall(VecDot(dX, dX, &dXdotdX)); 272 273 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 274 snes->norm = fnorm; 275 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 276 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 277 278 /* test convergence */ 279 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 280 PetscCall(SNESMonitor(snes, 0, fnorm)); 281 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 282 283 /* Call general purpose update function */ 284 PetscTryTypeMethod(snes, update, snes->iter); 285 286 /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 287 288 for (i = 1; i < maxits + 1; i++) { 289 /* some update types require the old update direction or conjugate direction */ 290 if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold)); 291 PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX)); 292 PetscCall(SNESLineSearchGetReason(linesearch, &lsresult)); 293 PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 294 if (lsresult) { 295 if (++snes->numFailures >= snes->maxFailures) { 296 snes->reason = SNES_DIVERGED_LINE_SEARCH; 297 PetscFunctionReturn(PETSC_SUCCESS); 298 } 299 } 300 if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 301 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 302 PetscFunctionReturn(PETSC_SUCCESS); 303 } 304 /* Monitor convergence */ 305 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 306 snes->iter = i; 307 snes->norm = fnorm; 308 snes->xnorm = xnorm; 309 snes->ynorm = ynorm; 310 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 311 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 312 313 /* Test for convergence */ 314 PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 315 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 316 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 317 318 /* Call general purpose update function */ 319 PetscTryTypeMethod(snes, update, snes->iter); 320 if (snes->npc) { 321 if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 322 PetscCall(SNESApplyNPC(snes, X, NULL, dX)); 323 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 324 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 325 snes->reason = SNES_DIVERGED_INNER; 326 PetscFunctionReturn(PETSC_SUCCESS); 327 } 328 PetscCall(VecCopy(dX, F)); 329 } else { 330 PetscCall(SNESApplyNPC(snes, X, F, dX)); 331 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 332 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 333 snes->reason = SNES_DIVERGED_INNER; 334 PetscFunctionReturn(PETSC_SUCCESS); 335 } 336 } 337 } else { 338 PetscCall(VecCopy(F, dX)); 339 } 340 341 /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 342 switch (ncg->type) { 343 case SNES_NCG_FR: /* Fletcher-Reeves */ 344 dXolddotdXold = dXdotdX; 345 PetscCall(VecDot(dX, dX, &dXdotdX)); 346 beta = PetscRealPart(dXdotdX / dXolddotdXold); 347 break; 348 case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 349 dXolddotdXold = dXdotdX; 350 PetscCall(VecDotBegin(dX, dX, &dXdotdXold)); 351 PetscCall(VecDotBegin(dXold, dX, &dXdotdXold)); 352 PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 353 PetscCall(VecDotEnd(dXold, dX, &dXdotdXold)); 354 beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 355 if (beta < 0.0) beta = 0.0; /* restart */ 356 break; 357 case SNES_NCG_HS: /* Hestenes-Stiefel */ 358 PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 359 PetscCall(VecDotBegin(dX, dXold, &dXdotdXold)); 360 PetscCall(VecDotBegin(lX, dX, &lXdotdX)); 361 PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 362 PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 363 PetscCall(VecDotEnd(dX, dXold, &dXdotdXold)); 364 PetscCall(VecDotEnd(lX, dX, &lXdotdX)); 365 PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 366 beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 367 break; 368 case SNES_NCG_DY: /* Dai-Yuan */ 369 PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 370 PetscCall(VecDotBegin(lX, dX, &lXdotdX)); 371 PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 372 PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 373 PetscCall(VecDotEnd(lX, dX, &lXdotdX)); 374 PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 375 beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX)); 376 break; 377 case SNES_NCG_CD: /* Conjugate Descent */ 378 PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 379 PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 380 PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 381 PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 382 beta = PetscRealPart(dXdotdX / lXdotdXold); 383 break; 384 } 385 if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta)); 386 PetscCall(VecAYPX(lX, beta, dX)); 387 } 388 PetscFunctionReturn(PETSC_SUCCESS); 389 } 390 391 /*MC 392 SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems. 393 394 Level: beginner 395 396 Options Database Keys: 397 + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp. 398 . -snes_linesearch_type <cp,l2,basic> - Line search type. 399 - -snes_ncg_monitor - Print the beta values used in the ncg iteration, . 400 401 Notes: 402 This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 403 gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 404 chooses the initial search direction as F(x) for the initial guess x. 405 406 Only supports left non-linear preconditioning. 407 408 Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then 409 `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR` 410 411 References: 412 . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 413 SIAM Review, 57(4), 2015 414 415 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()` 416 M*/ 417 PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 418 { 419 SNES_NCG *neP; 420 421 PetscFunctionBegin; 422 snes->ops->destroy = SNESDestroy_NCG; 423 snes->ops->setup = SNESSetUp_NCG; 424 snes->ops->setfromoptions = SNESSetFromOptions_NCG; 425 snes->ops->view = SNESView_NCG; 426 snes->ops->solve = SNESSolve_NCG; 427 snes->ops->reset = SNESReset_NCG; 428 429 snes->usesksp = PETSC_FALSE; 430 snes->usesnpc = PETSC_TRUE; 431 snes->npcside = PC_LEFT; 432 433 snes->alwayscomputesfinalresidual = PETSC_TRUE; 434 435 if (!snes->tolerancesset) { 436 snes->max_funcs = 30000; 437 snes->max_its = 10000; 438 snes->stol = 1e-20; 439 } 440 441 PetscCall(PetscNew(&neP)); 442 snes->data = (void *)neP; 443 neP->monitor = NULL; 444 neP->type = SNES_NCG_PRP; 445 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG)); 446 PetscFunctionReturn(PETSC_SUCCESS); 447 } 448