1 #include <private/snesimpl.h> 2 3 4 5 typedef struct { 6 Vec *dX; /* The change in X */ 7 Vec *dF; /* The change in F */ 8 PetscInt m; /* the number of kept previous steps */ 9 PetscScalar *alpha; 10 PetscScalar *beta; 11 PetscScalar *rho; 12 PetscViewer monitor; 13 PetscReal powell_gamma; /* Powell angle restart condition */ 14 PetscReal powell_downhill; /* Powell descent restart condition */ 15 PetscReal scaling; /* scaling of H0 */ 16 } SNES_QN; 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "LBGFSApplyJinv_Private" 20 PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec D, Vec Y) { 21 22 PetscErrorCode ierr; 23 24 SNES_QN *qn = (SNES_QN*)snes->data; 25 26 Vec *dX = qn->dX; 27 Vec *dF = qn->dF; 28 29 PetscScalar *alpha = qn->alpha; 30 PetscScalar *beta = qn->beta; 31 PetscScalar *rho = qn->rho; 32 33 PetscInt k, i; 34 PetscInt m = qn->m; 35 PetscScalar t; 36 PetscInt l = m; 37 38 PetscFunctionBegin; 39 40 ierr = VecCopy(D, Y);CHKERRQ(ierr); 41 42 if (it < m) l = it; 43 44 /* outward recursion starting at iteration k's update and working back */ 45 for (i = 0; i < l; i++) { 46 k = (it - i - 1) % l; 47 ierr = VecDot(dX[k], Y, &t);CHKERRQ(ierr); 48 alpha[k] = t*rho[k]; 49 if (qn->monitor) { 50 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 51 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 52 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 53 } 54 ierr = VecAXPY(Y, -alpha[k], dF[k]);CHKERRQ(ierr); 55 } 56 57 ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 58 59 /* inward recursion starting at the first update and working forward */ 60 for (i = 0; i < l; i++) { 61 k = (it + i - l) % l; 62 ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 63 beta[k] = rho[k]*t; 64 ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]); 65 if (qn->monitor) { 66 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 67 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 68 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 69 } 70 } 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "SNESSolve_QN" 76 static PetscErrorCode SNESSolve_QN(SNES snes) 77 { 78 79 PetscErrorCode ierr; 80 SNES_QN *qn = (SNES_QN*) snes->data; 81 82 Vec X, Xold; 83 Vec F, G, B; 84 Vec W, Y, FPC, Fold; 85 SNESConvergedReason reason; 86 PetscInt i, i_r, k; 87 88 PetscReal fnorm, xnorm = 0, ynorm, gnorm; 89 PetscInt m = qn->m; 90 PetscBool lssucceed; 91 92 PetscScalar rhosc; 93 94 Vec *dX = qn->dX; 95 Vec *dF = qn->dF; 96 PetscScalar *rho = qn->rho; 97 PetscScalar FolddotF, FolddotFold, FdotF, YdotF, a; 98 99 /* basically just a regular newton's method except for the application of the jacobian */ 100 PetscFunctionBegin; 101 102 X = snes->vec_sol; /* solution vector */ 103 F = snes->vec_func; /* residual vector */ 104 Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 105 B = snes->vec_rhs; 106 G = snes->work[0]; 107 W = snes->work[1]; 108 Xold = snes->work[2]; 109 110 /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 111 Fold = snes->work[3]; 112 113 snes->reason = SNES_CONVERGED_ITERATING; 114 115 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 116 snes->iter = 0; 117 snes->norm = 0.; 118 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 119 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 120 if (snes->domainerror) { 121 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 122 PetscFunctionReturn(0); 123 } 124 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 125 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 126 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 127 snes->norm = fnorm; 128 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 129 SNESLogConvHistory(snes,fnorm,0); 130 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 131 132 /* set parameter for default relative tolerance convergence test */ 133 snes->ttol = fnorm*snes->rtol; 134 /* test convergence */ 135 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 136 if (snes->reason) PetscFunctionReturn(0); 137 138 /* initialize the search direction as steepest descent */ 139 ierr = VecCopy(F, Y);CHKERRQ(ierr); 140 141 for(i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 142 /* line search for lambda */ 143 ynorm = 1; gnorm = fnorm; 144 ierr = VecCopy(F, Fold);CHKERRQ(ierr); 145 ierr = VecCopy(X, Xold);CHKERRQ(ierr); 146 ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 147 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 148 if (snes->domainerror) { 149 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 150 PetscFunctionReturn(0); 151 } 152 if (!lssucceed) { 153 if (++snes->numFailures >= snes->maxFailures) { 154 snes->reason = SNES_DIVERGED_LINE_SEARCH; 155 break; 156 } 157 } 158 159 /* Update function and solution vectors */ 160 fnorm = gnorm; 161 ierr = VecCopy(G,F);CHKERRQ(ierr); 162 ierr = VecCopy(W,X);CHKERRQ(ierr); 163 164 if (snes->pc) { 165 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 166 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 167 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 168 snes->reason = SNES_DIVERGED_INNER; 169 PetscFunctionReturn(0); 170 } 171 ierr = SNESGetFunction(snes->pc, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 172 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 173 ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr); 174 } 175 176 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 177 178 ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 179 ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 180 181 SNESLogConvHistory(snes,snes->norm,snes->iter); 182 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 183 /* set parameter for default relative tolerance convergence test */ 184 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 185 if (snes->reason) PetscFunctionReturn(0); 186 187 /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 188 ierr = VecDot(Fold, Fold, &FolddotFold);CHKERRQ(ierr); 189 ierr = VecDot(Fold, F, &FolddotF);CHKERRQ(ierr); 190 ierr = VecDot(F, F, &FdotF);CHKERRQ(ierr); 191 ierr = VecDot(Y, F, &YdotF);CHKERRQ(ierr); 192 if (PetscAbs(PetscRealPart(FolddotF)) > qn->powell_gamma*PetscAbs(PetscRealPart(FolddotFold))) { 193 if (qn->monitor) { 194 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e|\n", k, PetscRealPart(FolddotF), qn->powell_gamma, PetscRealPart(FolddotFold));CHKERRQ(ierr); 196 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 197 } 198 i_r = -1; 199 /* general purpose update */ 200 if (snes->ops->update) { 201 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 202 } 203 ierr = VecCopy(F, Y);CHKERRQ(ierr); 204 } else { 205 /* set the differences */ 206 k = i_r % m; 207 ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 208 ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 209 ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 210 ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 211 ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 212 213 /* set scaling to be shanno scaling */ 214 rho[k] = 1. / rhosc; 215 ierr = VecDot(dF[k], dF[k], &a);CHKERRQ(ierr); 216 qn->scaling = PetscRealPart(rhosc) / PetscRealPart(a); 217 218 /* general purpose update */ 219 if (snes->ops->update) { 220 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 221 } 222 /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 223 ierr = LBGFSApplyJinv_Private(snes, i_r+1, F, Y);CHKERRQ(ierr); 224 } 225 } 226 if (i == snes->max_its) { 227 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 228 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 229 } 230 PetscFunctionReturn(0); 231 } 232 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "SNESSetUp_QN" 236 static PetscErrorCode SNESSetUp_QN(SNES snes) 237 { 238 SNES_QN *qn = (SNES_QN*)snes->data; 239 PetscErrorCode ierr; 240 PetscFunctionBegin; 241 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 242 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 243 ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 244 ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "SNESReset_QN" 250 static PetscErrorCode SNESReset_QN(SNES snes) 251 { 252 PetscErrorCode ierr; 253 SNES_QN *qn; 254 PetscFunctionBegin; 255 if (snes->data) { 256 qn = (SNES_QN*)snes->data; 257 if (qn->dX) { 258 ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 259 } 260 if (qn->dF) { 261 ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 262 } 263 ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 264 } 265 PetscFunctionReturn(0); 266 } 267 268 #undef __FUNCT__ 269 #define __FUNCT__ "SNESDestroy_QN" 270 static PetscErrorCode SNESDestroy_QN(SNES snes) 271 { 272 PetscErrorCode ierr; 273 PetscFunctionBegin; 274 ierr = SNESReset_QN(snes);CHKERRQ(ierr); 275 ierr = PetscFree(snes->data);CHKERRQ(ierr); 276 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"","",PETSC_NULL);CHKERRQ(ierr); 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "SNESSetFromOptions_QN" 282 static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 283 { 284 285 PetscErrorCode ierr; 286 SNES_QN *qn; 287 PetscBool monflg = PETSC_FALSE; 288 PetscFunctionBegin; 289 290 qn = (SNES_QN*)snes->data; 291 292 ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 293 ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-BFGS methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 294 ierr = PetscOptionsReal("-snes_qn_powell_gamma", "Powell angle tolerance", "SNES", qn->powell_gamma, &qn->powell_gamma, PETSC_NULL);CHKERRQ(ierr); 295 ierr = PetscOptionsReal("-snes_qn_powell_downhill", "Powell descent tolerance", "SNES", qn->powell_downhill, &qn->powell_downhill, PETSC_NULL);CHKERRQ(ierr); 296 ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNES", monflg, &monflg, PETSC_NULL);CHKERRQ(ierr); 297 ierr = PetscOptionsTail();CHKERRQ(ierr); 298 if (monflg) { 299 qn->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 300 } 301 PetscFunctionReturn(0); 302 } 303 304 EXTERN_C_BEGIN 305 #undef __FUNCT__ 306 #define __FUNCT__ "SNESLineSearchSetType_QN" 307 PetscErrorCode SNESLineSearchSetType_QN(SNES snes, SNESLineSearchType type) 308 { 309 PetscErrorCode ierr; 310 PetscFunctionBegin; 311 312 switch (type) { 313 case SNES_LS_BASIC: 314 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 315 break; 316 case SNES_LS_BASIC_NONORMS: 317 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 318 break; 319 case SNES_LS_QUADRATIC: 320 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 321 break; 322 case SNES_LS_SECANT: 323 ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 324 break; 325 default: 326 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 327 break; 328 } 329 snes->ls_type = type; 330 PetscFunctionReturn(0); 331 } 332 EXTERN_C_END 333 334 335 /* -------------------------------------------------------------------------- */ 336 /*MC 337 SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 338 339 Options Database: 340 341 + -snes_qn_m - Number of past states saved for the L-Broyden methods. 342 - -snes_qn_monitor - Monitors the quasi-newton jacobian. 343 344 Notes: This implements the L-BFGS algorithm for the solution of F(x) = b using previous change in F(x) and x to 345 form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 346 generalized to implement several limited-memory Broyden methods. 347 348 References: 349 350 L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 351 International Journal of Computer Mathematics, vol. 86, 2009. 352 353 354 Level: beginner 355 356 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 357 358 M*/ 359 EXTERN_C_BEGIN 360 #undef __FUNCT__ 361 #define __FUNCT__ "SNESCreate_QN" 362 PetscErrorCode SNESCreate_QN(SNES snes) 363 { 364 365 PetscErrorCode ierr; 366 SNES_QN *qn; 367 368 PetscFunctionBegin; 369 snes->ops->setup = SNESSetUp_QN; 370 snes->ops->solve = SNESSolve_QN; 371 snes->ops->destroy = SNESDestroy_QN; 372 snes->ops->setfromoptions = SNESSetFromOptions_QN; 373 snes->ops->view = 0; 374 snes->ops->reset = SNESReset_QN; 375 376 snes->usespc = PETSC_TRUE; 377 snes->usesksp = PETSC_FALSE; 378 379 snes->max_funcs = 30000; 380 snes->max_its = 10000; 381 382 ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 383 snes->data = (void *) qn; 384 qn->m = 10; 385 qn->scaling = 1.0; 386 qn->dX = PETSC_NULL; 387 qn->dF = PETSC_NULL; 388 qn->monitor = PETSC_NULL; 389 qn->powell_gamma = 0.8; 390 qn->powell_downhill = 0.4; 391 392 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_QN",SNESLineSearchSetType_QN);CHKERRQ(ierr); 393 ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr); 394 395 PetscFunctionReturn(0); 396 } 397 EXTERN_C_END 398