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