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