1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdm.h> 3 4 #define H(i, j) qn->dXdFmat[i * qn->m + j] 5 6 const char *const SNESQNScaleTypes[] = {"DEFAULT", "NONE", "SCALAR", "DIAGONAL", "JACOBIAN", "SNESQNScaleType", "SNES_QN_SCALING_", NULL}; 7 const char *const SNESQNRestartTypes[] = {"DEFAULT", "NONE", "POWELL", "PERIODIC", "SNESQNRestartType", "SNES_QN_RESTART_", NULL}; 8 const char *const SNESQNTypes[] = {"LBFGS", "BROYDEN", "BADBROYDEN", "SNESQNType", "SNES_QN_", NULL}; 9 10 typedef struct { 11 Mat B; /* Quasi-Newton approximation Matrix (MATLMVM) */ 12 PetscInt m; /* The number of kept previous steps */ 13 PetscReal *lambda; /* The line search history of the method */ 14 PetscBool monflg; 15 PetscViewer monitor; 16 PetscReal powell_gamma; /* Powell angle restart condition */ 17 PetscReal scaling; /* scaling of H0 */ 18 SNESQNType type; /* the type of quasi-newton method used */ 19 SNESQNScaleType scale_type; /* the type of scaling used */ 20 SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 21 } SNES_QN; 22 23 static PetscErrorCode SNESSolve_QN(SNES snes) 24 { 25 SNES_QN *qn = (SNES_QN *)snes->data; 26 Vec X, Xold; 27 Vec F, W; 28 Vec Y, D, Dold; 29 PetscInt i, i_r; 30 PetscReal fnorm, xnorm, ynorm, gnorm; 31 SNESLineSearchReason lssucceed; 32 PetscBool badstep, powell, periodic; 33 PetscScalar DolddotD, DolddotDold; 34 SNESConvergedReason reason; 35 36 /* basically just a regular newton's method except for the application of the Jacobian */ 37 38 PetscFunctionBegin; 39 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); 40 41 PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 42 F = snes->vec_func; /* residual vector */ 43 Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 44 W = snes->work[3]; 45 X = snes->vec_sol; /* solution vector */ 46 Xold = snes->work[0]; 47 48 /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 49 D = snes->work[1]; 50 Dold = snes->work[2]; 51 52 snes->reason = SNES_CONVERGED_ITERATING; 53 54 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 55 snes->iter = 0; 56 snes->norm = 0.; 57 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 58 59 if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 60 PetscCall(SNESApplyNPC(snes, X, NULL, F)); 61 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 62 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 63 snes->reason = SNES_DIVERGED_INNER; 64 PetscFunctionReturn(PETSC_SUCCESS); 65 } 66 PetscCall(VecNorm(F, NORM_2, &fnorm)); 67 } else { 68 if (!snes->vec_func_init_set) { 69 PetscCall(SNESComputeFunction(snes, X, F)); 70 } else snes->vec_func_init_set = PETSC_FALSE; 71 72 PetscCall(VecNorm(F, NORM_2, &fnorm)); 73 SNESCheckFunctionNorm(snes, fnorm); 74 } 75 if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 76 PetscCall(SNESApplyNPC(snes, X, F, D)); 77 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 78 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 79 snes->reason = SNES_DIVERGED_INNER; 80 PetscFunctionReturn(PETSC_SUCCESS); 81 } 82 } else { 83 PetscCall(VecCopy(F, D)); 84 } 85 86 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 87 snes->norm = fnorm; 88 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 89 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 90 PetscCall(SNESMonitor(snes, 0, fnorm)); 91 92 /* test convergence */ 93 PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 94 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 95 96 if (snes->npc && snes->npcside == PC_RIGHT) { 97 PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0)); 98 PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 99 PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0)); 100 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 101 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 102 snes->reason = SNES_DIVERGED_INNER; 103 PetscFunctionReturn(PETSC_SUCCESS); 104 } 105 PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 106 PetscCall(VecCopy(F, D)); 107 } 108 109 /* general purpose update */ 110 PetscTryTypeMethod(snes, update, snes->iter); 111 112 /* scale the initial update */ 113 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 114 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 115 SNESCheckJacobianDomainerror(snes); 116 PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 117 PetscCall(MatLMVMSetJ0KSP(qn->B, snes->ksp)); 118 } 119 120 for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 121 /* update QN approx and calculate step */ 122 PetscCall(MatLMVMUpdate(qn->B, X, D)); 123 PetscCall(MatSolve(qn->B, D, Y)); 124 125 /* line search for lambda */ 126 ynorm = 1; 127 gnorm = fnorm; 128 PetscCall(VecCopy(D, Dold)); 129 PetscCall(VecCopy(X, Xold)); 130 PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y)); 131 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 132 PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 133 PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm)); 134 badstep = PETSC_FALSE; 135 if (lssucceed) { 136 if (++snes->numFailures >= snes->maxFailures) { 137 snes->reason = SNES_DIVERGED_LINE_SEARCH; 138 break; 139 } 140 badstep = PETSC_TRUE; 141 } 142 143 /* convergence monitoring */ 144 PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)fnorm, (double)gnorm, (double)ynorm, (int)lssucceed)); 145 146 if (snes->npc && snes->npcside == PC_RIGHT) { 147 PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, 0, 0)); 148 PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 149 PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, 0, 0)); 150 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 151 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 152 snes->reason = SNES_DIVERGED_INNER; 153 PetscFunctionReturn(PETSC_SUCCESS); 154 } 155 PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 156 } 157 158 PetscCall(SNESSetIterationNumber(snes, i + 1)); 159 snes->norm = fnorm; 160 snes->xnorm = xnorm; 161 snes->ynorm = ynorm; 162 163 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter)); 164 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 165 166 /* set parameter for default relative tolerance convergence test */ 167 PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 168 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 169 if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 170 PetscCall(SNESApplyNPC(snes, X, F, D)); 171 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 172 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 173 snes->reason = SNES_DIVERGED_INNER; 174 PetscFunctionReturn(PETSC_SUCCESS); 175 } 176 } else { 177 PetscCall(VecCopy(F, D)); 178 } 179 180 /* general purpose update */ 181 PetscTryTypeMethod(snes, update, snes->iter); 182 183 /* restart conditions */ 184 powell = PETSC_FALSE; 185 if (qn->restart_type == SNES_QN_RESTART_POWELL && i_r > 1) { 186 /* check restart by Powell's Criterion: |F^T H_0 Fold| > powell_gamma * |Fold^T H_0 Fold| */ 187 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 188 PetscCall(MatMult(snes->jacobian, Dold, W)); 189 } else { 190 PetscCall(VecCopy(Dold, W)); 191 } 192 PetscCall(VecDotBegin(W, Dold, &DolddotDold)); 193 PetscCall(VecDotBegin(W, D, &DolddotD)); 194 PetscCall(VecDotEnd(W, Dold, &DolddotDold)); 195 PetscCall(VecDotEnd(W, D, &DolddotD)); 196 if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma * PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 197 } 198 periodic = PETSC_FALSE; 199 if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 200 if (i_r > qn->m - 1) periodic = PETSC_TRUE; 201 } 202 /* restart if either powell or periodic restart is satisfied. */ 203 if (badstep || powell || periodic) { 204 if (qn->monflg) { 205 PetscCall(PetscViewerASCIIAddTab(qn->monitor, ((PetscObject)snes)->tablevel + 2)); 206 if (powell) { 207 PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Powell restart! |%14.12e| > %6.4f*|%14.12e| i_r = %" PetscInt_FMT "\n", (double)PetscRealPart(DolddotD), (double)qn->powell_gamma, (double)PetscRealPart(DolddotDold), i_r)); 208 } else { 209 PetscCall(PetscViewerASCIIPrintf(qn->monitor, "Periodic restart! i_r = %" PetscInt_FMT "\n", i_r)); 210 } 211 PetscCall(PetscViewerASCIISubtractTab(qn->monitor, ((PetscObject)snes)->tablevel + 2)); 212 } 213 i_r = -1; 214 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 215 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 216 SNESCheckJacobianDomainerror(snes); 217 } 218 PetscCall(MatLMVMReset(qn->B, PETSC_FALSE)); 219 } 220 } 221 if (i == snes->max_its) { 222 PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its)); 223 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 224 } 225 PetscFunctionReturn(PETSC_SUCCESS); 226 } 227 228 static PetscErrorCode SNESSetUp_QN(SNES snes) 229 { 230 SNES_QN *qn = (SNES_QN *)snes->data; 231 DM dm; 232 PetscInt n, N; 233 234 PetscFunctionBegin; 235 236 if (!snes->vec_sol) { 237 PetscCall(SNESGetDM(snes, &dm)); 238 PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol)); 239 } 240 PetscCall(SNESSetWorkVecs(snes, 4)); 241 242 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) PetscCall(SNESSetUpMatrices(snes)); 243 if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 244 245 /* set method defaults */ 246 if (qn->scale_type == SNES_QN_SCALE_DEFAULT) { 247 if (qn->type == SNES_QN_BADBROYDEN) { 248 qn->scale_type = SNES_QN_SCALE_NONE; 249 } else { 250 qn->scale_type = SNES_QN_SCALE_SCALAR; 251 } 252 } 253 if (qn->restart_type == SNES_QN_RESTART_DEFAULT) { 254 if (qn->type == SNES_QN_LBFGS) { 255 qn->restart_type = SNES_QN_RESTART_POWELL; 256 } else { 257 qn->restart_type = SNES_QN_RESTART_PERIODIC; 258 } 259 } 260 261 /* Set up the LMVM matrix */ 262 switch (qn->type) { 263 case SNES_QN_BROYDEN: 264 PetscCall(MatSetType(qn->B, MATLMVMBROYDEN)); 265 qn->scale_type = SNES_QN_SCALE_NONE; 266 break; 267 case SNES_QN_BADBROYDEN: 268 PetscCall(MatSetType(qn->B, MATLMVMBADBROYDEN)); 269 qn->scale_type = SNES_QN_SCALE_NONE; 270 break; 271 default: 272 PetscCall(MatSetType(qn->B, MATLMVMBFGS)); 273 switch (qn->scale_type) { 274 case SNES_QN_SCALE_NONE: 275 PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_NONE)); 276 break; 277 case SNES_QN_SCALE_SCALAR: 278 PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_SCALAR)); 279 break; 280 case SNES_QN_SCALE_JACOBIAN: 281 PetscCall(MatLMVMSymBroydenSetScaleType(qn->B, MAT_LMVM_SYMBROYDEN_SCALE_USER)); 282 break; 283 case SNES_QN_SCALE_DIAGONAL: 284 case SNES_QN_SCALE_DEFAULT: 285 default: 286 break; 287 } 288 break; 289 } 290 PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 291 PetscCall(VecGetSize(snes->vec_sol, &N)); 292 PetscCall(MatSetSizes(qn->B, n, n, N, N)); 293 PetscCall(MatSetUp(qn->B)); 294 PetscCall(MatLMVMReset(qn->B, PETSC_TRUE)); 295 PetscCall(MatLMVMSetHistorySize(qn->B, qn->m)); 296 PetscCall(MatLMVMAllocate(qn->B, snes->vec_sol, snes->vec_func)); 297 PetscFunctionReturn(PETSC_SUCCESS); 298 } 299 300 static PetscErrorCode SNESReset_QN(SNES snes) 301 { 302 SNES_QN *qn; 303 304 PetscFunctionBegin; 305 if (snes->data) { 306 qn = (SNES_QN *)snes->data; 307 PetscCall(MatDestroy(&qn->B)); 308 } 309 PetscFunctionReturn(PETSC_SUCCESS); 310 } 311 312 static PetscErrorCode SNESDestroy_QN(SNES snes) 313 { 314 PetscFunctionBegin; 315 PetscCall(SNESReset_QN(snes)); 316 PetscCall(PetscFree(snes->data)); 317 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", NULL)); 318 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", NULL)); 319 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", NULL)); 320 PetscFunctionReturn(PETSC_SUCCESS); 321 } 322 323 static PetscErrorCode SNESSetFromOptions_QN(SNES snes, PetscOptionItems *PetscOptionsObject) 324 { 325 SNES_QN *qn = (SNES_QN *)snes->data; 326 PetscBool flg; 327 SNESLineSearch linesearch; 328 SNESQNRestartType rtype = qn->restart_type; 329 SNESQNScaleType stype = qn->scale_type; 330 SNESQNType qtype = qn->type; 331 332 PetscFunctionBegin; 333 PetscOptionsHeadBegin(PetscOptionsObject, "SNES QN options"); 334 PetscCall(PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNESQN", qn->m, &qn->m, NULL)); 335 PetscCall(PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL)); 336 PetscCall(PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", qn->monflg, &qn->monflg, NULL)); 337 PetscCall(PetscOptionsEnum("-snes_qn_scale_type", "Scaling type", "SNESQNSetScaleType", SNESQNScaleTypes, (PetscEnum)stype, (PetscEnum *)&stype, &flg)); 338 if (flg) PetscCall(SNESQNSetScaleType(snes, stype)); 339 340 PetscCall(PetscOptionsEnum("-snes_qn_restart_type", "Restart type", "SNESQNSetRestartType", SNESQNRestartTypes, (PetscEnum)rtype, (PetscEnum *)&rtype, &flg)); 341 if (flg) PetscCall(SNESQNSetRestartType(snes, rtype)); 342 343 PetscCall(PetscOptionsEnum("-snes_qn_type", "Quasi-Newton update type", "", SNESQNTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, &flg)); 344 if (flg) PetscCall(SNESQNSetType(snes, qtype)); 345 PetscCall(MatSetFromOptions(qn->B)); 346 PetscOptionsHeadEnd(); 347 if (!snes->linesearch) { 348 PetscCall(SNESGetLineSearch(snes, &linesearch)); 349 if (!((PetscObject)linesearch)->type_name) { 350 if (qn->type == SNES_QN_LBFGS) { 351 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP)); 352 } else if (qn->type == SNES_QN_BROYDEN) { 353 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 354 } else { 355 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 356 } 357 } 358 } 359 if (qn->monflg) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes), &qn->monitor)); 360 PetscFunctionReturn(PETSC_SUCCESS); 361 } 362 363 static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer) 364 { 365 SNES_QN *qn = (SNES_QN *)snes->data; 366 PetscBool iascii; 367 368 PetscFunctionBegin; 369 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 370 if (iascii) { 371 PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, restart type is %s, scale type is %s\n", SNESQNTypes[qn->type], SNESQNRestartTypes[qn->restart_type], SNESQNScaleTypes[qn->scale_type])); 372 PetscCall(PetscViewerASCIIPrintf(viewer, " Stored subspace size: %" PetscInt_FMT "\n", qn->m)); 373 } 374 PetscFunctionReturn(PETSC_SUCCESS); 375 } 376 377 /*@ 378 SNESQNSetRestartType - Sets the restart type for `SNESQN`. 379 380 Logically Collective 381 382 Input Parameters: 383 + snes - the iterative context 384 - rtype - restart type 385 386 Options Database Keys: 387 + -snes_qn_restart_type <powell,periodic,none> - set the restart type 388 - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic 389 390 Level: intermediate 391 392 `SNESQNRestartType`s: 393 + `SNES_QN_RESTART_NONE` - never restart 394 . `SNES_QN_RESTART_POWELL` - restart based upon descent criteria 395 - `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations 396 397 .seealso: `SNESQN`, `SNESQNRestartType`, `SNES_QN_RESTART_NONE`, `SNES_QN_RESTART_POWELL`, `SNES_QN_RESTART_PERIODIC` 398 @*/ 399 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 400 { 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 403 PetscTryMethod(snes, "SNESQNSetRestartType_C", (SNES, SNESQNRestartType), (snes, rtype)); 404 PetscFunctionReturn(PETSC_SUCCESS); 405 } 406 407 /*@ 408 SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in `SNESQN`. 409 410 Logically Collective 411 412 Input Parameters: 413 + snes - the nonlinear solver context 414 - stype - scale type 415 416 Options Database Key: 417 . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type 418 419 Level: intermediate 420 421 `SNESQNScaleType`s: 422 + `SNES_QN_SCALE_NONE` - don't scale the problem 423 . `SNES_QN_SCALE_SCALAR` - use Shanno scaling 424 . `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available 425 - `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration 426 of QN and at ever restart. 427 428 .seealso: `SNES`, `SNESQN`, `SNESLineSearch`, `SNESQNScaleType`, `SNESSetJacobian()` 429 @*/ 430 431 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 432 { 433 PetscFunctionBegin; 434 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 435 PetscTryMethod(snes, "SNESQNSetScaleType_C", (SNES, SNESQNScaleType), (snes, stype)); 436 PetscFunctionReturn(PETSC_SUCCESS); 437 } 438 439 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 440 { 441 SNES_QN *qn = (SNES_QN *)snes->data; 442 443 PetscFunctionBegin; 444 qn->scale_type = stype; 445 if (stype == SNES_QN_SCALE_JACOBIAN) snes->usesksp = PETSC_TRUE; 446 PetscFunctionReturn(PETSC_SUCCESS); 447 } 448 449 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 450 { 451 SNES_QN *qn = (SNES_QN *)snes->data; 452 453 PetscFunctionBegin; 454 qn->restart_type = rtype; 455 PetscFunctionReturn(PETSC_SUCCESS); 456 } 457 458 /*@ 459 SNESQNSetType - Sets the quasi-Newton variant to be used in `SNESQN`. 460 461 Logically Collective 462 463 Input Parameters: 464 + snes - the iterative context 465 - qtype - variant type 466 467 Options Database Key: 468 . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type 469 470 Level: beginner 471 472 `SNESQNType`s: 473 + `SNES_QN_LBFGS` - LBFGS variant 474 . `SNES_QN_BROYDEN` - Broyden variant 475 - `SNES_QN_BADBROYDEN` - Bad Broyden variant 476 477 .seealso: `SNESQN`, `SNES_QN_LBFGS`, `SNES_QN_BROYDEN`, `SNES_QN_BADBROYDEN`, `SNESQNType`, `TAOLMVM`, `TAOBLMVM` 478 @*/ 479 480 PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype) 481 { 482 PetscFunctionBegin; 483 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 484 PetscTryMethod(snes, "SNESQNSetType_C", (SNES, SNESQNType), (snes, qtype)); 485 PetscFunctionReturn(PETSC_SUCCESS); 486 } 487 488 PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype) 489 { 490 SNES_QN *qn = (SNES_QN *)snes->data; 491 492 PetscFunctionBegin; 493 qn->type = qtype; 494 PetscFunctionReturn(PETSC_SUCCESS); 495 } 496 497 /*MC 498 SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 499 500 Options Database Keys: 501 + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods. 502 . -snes_qn_restart_type <powell,periodic,none> - set the restart type 503 . -snes_qn_powell_gamma - Angle condition for restart. 504 . -snes_qn_powell_descent - Descent condition for restart. 505 . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type 506 . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - scaling performed on inner Jacobian 507 . -snes_linesearch_type <cp, l2, basic> - Type of line search. 508 - -snes_qn_monitor - Monitors the quasi-newton Jacobian. 509 510 References: 511 + * - Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995. 512 . * - R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods, 513 Technical Report, Northwestern University, June 1992. 514 . * - Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE 515 Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 516 . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 517 SIAM Review, 57(4), 2015 518 . * - Griewank, Andreas. "Broyden updating, the good and the bad!." Doc. Math (2012): 301-315. 519 . * - Gilbert, Jean Charles, and Claude Lemar{\'e}chal. "Some numerical experiments with variable-storage quasi-Newton algorithms." 520 Mathematical programming 45.1-3 (1989): 407-435. 521 - * - Dener A., Munson T. "Accelerating Limited-Memory Quasi-Newton Convergence for Large-Scale Optimization" 522 Computational Science - ICCS 2019. ICCS 2019. Lecture Notes in Computer Science, vol 11538. Springer, Cham 523 524 Level: beginner 525 526 Notes: 527 This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 528 previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 529 updates. 530 531 When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 532 these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 533 iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed, 534 perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 535 536 Uses left nonlinear preconditioning by default. 537 538 .seealso: `SNESQNRestartType`, `SNESQNSetRestartType()`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, 539 `SNESQNScaleType`, `SNESQNSetScaleType()`, `SNESQNSetType`, `SNESQNSetType ()` 540 M*/ 541 PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 542 { 543 SNES_QN *qn; 544 const char *optionsprefix; 545 546 PetscFunctionBegin; 547 snes->ops->setup = SNESSetUp_QN; 548 snes->ops->solve = SNESSolve_QN; 549 snes->ops->destroy = SNESDestroy_QN; 550 snes->ops->setfromoptions = SNESSetFromOptions_QN; 551 snes->ops->view = SNESView_QN; 552 snes->ops->reset = SNESReset_QN; 553 554 snes->npcside = PC_LEFT; 555 556 snes->usesnpc = PETSC_TRUE; 557 snes->usesksp = PETSC_FALSE; 558 559 snes->alwayscomputesfinalresidual = PETSC_TRUE; 560 561 if (!snes->tolerancesset) { 562 snes->max_funcs = 30000; 563 snes->max_its = 10000; 564 } 565 566 PetscCall(PetscNew(&qn)); 567 snes->data = (void *)qn; 568 qn->m = 10; 569 qn->scaling = 1.0; 570 qn->monitor = NULL; 571 qn->monflg = PETSC_FALSE; 572 qn->powell_gamma = 0.9999; 573 qn->scale_type = SNES_QN_SCALE_DEFAULT; 574 qn->restart_type = SNES_QN_RESTART_DEFAULT; 575 qn->type = SNES_QN_LBFGS; 576 577 PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &qn->B)); 578 PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 579 PetscCall(MatSetOptionsPrefix(qn->B, optionsprefix)); 580 581 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetScaleType_C", SNESQNSetScaleType_QN)); 582 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetRestartType_C", SNESQNSetRestartType_QN)); 583 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESQNSetType_C", SNESQNSetType_QN)); 584 PetscFunctionReturn(PETSC_SUCCESS); 585 } 586