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