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