1 #include <private/snesimpl.h> 2 3 typedef struct { 4 Vec * dX; /* The change in X */ 5 Vec * dF; /* The change in F */ 6 PetscInt m; /* the number of kept previous steps */ 7 PetscScalar * alpha; 8 PetscScalar * beta; 9 PetscScalar * rho; 10 } QNContext; 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "LBGFSApplyJinv_Private" 14 PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) { 15 16 PetscErrorCode ierr; 17 18 QNContext * qn = (QNContext *)snes->data; 19 20 Vec * dX = qn->dX; 21 Vec * dF = qn->dF; 22 23 PetscScalar * alpha = qn->alpha; 24 PetscScalar * beta = qn->beta; 25 PetscScalar * rho = qn->rho; 26 27 PetscInt k, i; 28 PetscInt m = qn->m; 29 PetscScalar t; 30 PetscInt l = m; 31 32 PetscFunctionBegin; 33 34 if (it < m) l = it; 35 36 ierr = VecCopy(g, z);CHKERRQ(ierr); 37 38 /* outward recursion starting at iteration k's update and working back */ 39 for (i = 0; i < l; i++) { 40 k = (it - i - 1) % l; 41 ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr); 42 alpha[k] = t*rho[k]; 43 ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha %g\n", k, alpha[k]);CHKERRQ(ierr); 44 ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr); 45 } 46 47 /* inner application of the initial inverse jacobian approximation */ 48 /* right now it's just the identity. Nothing needs to go here. */ 49 50 /* inward recursion starting at the first update and working forward*/ 51 for (i = 0; i < l; i++) { 52 k = (it + i - l) % l; 53 ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr); 54 beta[k] = rho[k]*t; 55 ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]); 56 ierr = PetscPrintf(PETSC_COMM_WORLD, "%d, alpha - beta %g\n", k, alpha[k] - beta[k]);CHKERRQ(ierr); 57 } 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "SNESLineSearchQuadratic_QN" 63 PetscErrorCode SNESLineSearchQuadratic_QN(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal xnorm,Vec G,Vec W,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 64 { 65 PetscInt i; 66 PetscReal alphas[3] = {0.0, 0.5, 1.0}; 67 PetscReal norms[3]; 68 PetscReal alpha,a,b; 69 PetscErrorCode ierr; 70 71 PetscFunctionBegin; 72 norms[0] = fnorm; 73 for(i=1; i < 3; ++i) { 74 ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 75 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 76 ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr); 77 } 78 for(i = 0; i < 3; ++i) { 79 norms[i] = PetscSqr(norms[i]); 80 } 81 /* Fit a quadratic: 82 If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 83 a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1) 84 b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1) 85 c = y_0 86 x_min = -b/2a 87 88 If we let x_{0,1,2} = 0, 0.5, 1.0 89 a = 2 y_2 - 4 y_1 + 2 y_0 90 b = -y_2 + 4 y_1 - 3 y_0 91 c = y_0 92 */ 93 a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1])); 94 b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]); 95 /* Check for positive a (concave up) */ 96 if (a >= 0.0) { 97 alpha = -b/(2.0*a); 98 alpha = PetscMin(alpha, alphas[2]); 99 alpha = PetscMax(alpha, alphas[0]); 100 } else { 101 alpha = 1.0; 102 } 103 if (snes->ls_monitor) { 104 ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 105 ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr); 106 ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 107 } 108 ierr = VecCopy(X, W);CHKERRQ(ierr); 109 ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr); 110 if (alpha != 1.0) { 111 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 112 ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 113 } else { 114 *gnorm = PetscSqrtReal(norms[2]); 115 } 116 if (alpha == 0.0) *flag = PETSC_FALSE; 117 else *flag = PETSC_TRUE; 118 PetscFunctionReturn(0); 119 } 120 121 #undef __FUNCT__ 122 #define __FUNCT__ "SNESSolve_QN" 123 static PetscErrorCode SNESSolve_QN(SNES snes) 124 { 125 126 PetscErrorCode ierr; 127 QNContext * qn = (QNContext*) snes->data; 128 129 Vec X, Xold; 130 Vec F, Fold, G; 131 Vec W, Y; 132 133 PetscInt i, k; 134 135 PetscReal fnorm, xnorm = 0, ynorm, gnorm; 136 PetscInt m = qn->m; 137 PetscBool lssucceed; 138 139 PetscScalar rhosc; 140 141 Vec * dX = qn->dX; 142 Vec * dF = qn->dF; 143 PetscScalar * rho = qn->rho; 144 145 /* basically just a regular newton's method except for the application of the jacobian */ 146 PetscFunctionBegin; 147 148 X = snes->vec_sol; /* solution vector */ 149 F = snes->vec_func; /* residual vector */ 150 Y = snes->vec_sol_update; /* search direction */ 151 G = snes->work[0]; 152 W = snes->work[1]; 153 Xold = snes->work[2]; 154 Fold = snes->work[3]; 155 156 snes->reason = SNES_CONVERGED_ITERATING; 157 158 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 159 snes->iter = 0; 160 snes->norm = 0.; 161 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 162 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 163 if (snes->domainerror) { 164 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 165 PetscFunctionReturn(0); 166 } 167 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 168 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 169 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 170 snes->norm = fnorm; 171 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 172 SNESLogConvHistory(snes,fnorm,0); 173 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 174 175 /* set parameter for default relative tolerance convergence test */ 176 snes->ttol = fnorm*snes->rtol; 177 /* test convergence */ 178 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 179 if (snes->reason) PetscFunctionReturn(0); 180 181 /* initialize the search direction as steepest descent */ 182 ierr = VecCopy(F, Y);CHKERRQ(ierr); 183 for(i = 0; i < snes->max_its; i++) { 184 /* line search for lambda */ 185 ynorm = 1; gnorm = fnorm; 186 ierr = VecCopy(F, Fold);CHKERRQ(ierr); 187 ierr = VecCopy(X, Xold);CHKERRQ(ierr); 188 ierr = VecScale(Y, -1.0);CHKERRQ(ierr); 189 ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 190 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); 191 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 192 if (snes->domainerror) { 193 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 194 PetscFunctionReturn(0); 195 } 196 if (!lssucceed) { 197 if (++snes->numFailures >= snes->maxFailures) { 198 snes->reason = SNES_DIVERGED_LINE_SEARCH; 199 break; 200 } 201 } 202 /* Update function and solution vectors */ 203 fnorm = gnorm; 204 ierr = VecCopy(G,F);CHKERRQ(ierr); 205 ierr = VecCopy(W,X);CHKERRQ(ierr); 206 207 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 208 snes->iter = i+1; 209 snes->norm = fnorm; 210 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 211 SNESLogConvHistory(snes,snes->norm,snes->iter); 212 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 213 /* set parameter for default relative tolerance convergence test */ 214 ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 215 if (snes->reason) PetscFunctionReturn(0); 216 217 /* set the differences */ 218 k = i % m; 219 ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 220 ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 221 ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 222 ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 223 ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 224 rho[k] = 1. / rhosc; 225 226 /* general purpose update */ 227 if (snes->ops->update) { 228 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 229 } 230 231 /* apply the current iteration of the approximate jacobian in order to get the next search direction*/ 232 ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr); 233 } 234 if (i == snes->max_its) { 235 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 236 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 237 } 238 PetscFunctionReturn(0); 239 } 240 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "SNESSetUp_QN" 244 static PetscErrorCode SNESSetUp_QN(SNES snes) 245 { 246 QNContext * qn = (QNContext *)snes->data; 247 PetscErrorCode ierr; 248 PetscFunctionBegin; 249 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 250 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 251 ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 252 ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 253 PetscFunctionReturn(0); 254 } 255 256 #undef __FUNCT__ 257 #define __FUNCT__ "SNESReset_QN" 258 static PetscErrorCode SNESReset_QN(SNES snes) 259 { 260 PetscErrorCode ierr; 261 QNContext * qn; 262 PetscFunctionBegin; 263 if (snes->data) { 264 qn = (QNContext *)snes->data; 265 if (qn->dX) { 266 ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 267 } 268 if (qn->dF) { 269 ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 270 } 271 ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 272 } 273 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 274 PetscFunctionReturn(0); 275 } 276 277 #undef __FUNCT__ 278 #define __FUNCT__ "SNESDestroy_QN" 279 static PetscErrorCode SNESDestroy_QN(SNES snes) 280 { 281 PetscErrorCode ierr; 282 PetscFunctionBegin; 283 ierr = SNESReset_QN(snes);CHKERRQ(ierr); 284 ierr = PetscFree(snes->data);CHKERRQ(ierr); 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "SNESSetFromOptions_QN" 290 static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 291 { 292 293 PetscErrorCode ierr; 294 QNContext * qn; 295 296 PetscFunctionBegin; 297 298 qn = (QNContext *)snes->data; 299 300 ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 301 ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 302 ierr = PetscOptionsTail();CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 /* -------------------------------------------------------------------------- */ 307 /*MC 308 SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 309 310 Options Database: 311 312 + -snes_qn_m - Number of past states saved for the L-Broyden methods. 313 314 Notes: This implements the L-BFGS algorithm for the solution of F(x) = 0 using previous change in F(x) and x to 315 form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 316 generalized to implement several limited-memory Broyden methods. 317 318 References: 319 320 L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 321 International Journal of Computer Mathematics, vol. 86, 2009. 322 323 324 Level: beginner 325 326 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 327 328 M*/ 329 EXTERN_C_BEGIN 330 #undef __FUNCT__ 331 #define __FUNCT__ "SNESCreate_QN" 332 PetscErrorCode SNESCreate_QN(SNES snes) 333 { 334 335 PetscErrorCode ierr; 336 QNContext * qn; 337 338 PetscFunctionBegin; 339 snes->ops->setup = SNESSetUp_QN; 340 snes->ops->solve = SNESSolve_QN; 341 snes->ops->destroy = SNESDestroy_QN; 342 snes->ops->setfromoptions = SNESSetFromOptions_QN; 343 snes->ops->view = 0; 344 snes->ops->reset = SNESReset_QN; 345 346 snes->usespc = PETSC_TRUE; 347 snes->usesksp = PETSC_FALSE; 348 349 ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 350 snes->data = (void *) qn; 351 qn->m = 10; 352 qn->dX = PETSC_NULL; 353 qn->dF = PETSC_NULL; 354 355 snes->ops->linesearchcubic = SNESLineSearchCubic; 356 snes->ops->linesearchquadratic = SNESLineSearchQuadratic_QN; 357 snes->ops->linesearchno = SNESLineSearchNo; 358 snes->ops->linesearchnonorms = SNESLineSearchNoNorms; 359 360 ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 361 362 PetscFunctionReturn(0); 363 } 364 EXTERN_C_END 365