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) % m; 41 /* k = (it + i - l) % m; */ 42 ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr); 43 alpha[k] = t*rho[k]; 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 - 1) % m; */ 53 k = (it + i - l) % m; 54 ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr); 55 beta[k] = rho[k]*t; 56 ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]); 57 } 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "SNESSolve_QN" 63 static PetscErrorCode SNESSolve_QN(SNES snes) 64 { 65 66 PetscErrorCode ierr; 67 QNContext * qn = (QNContext*) snes->data; 68 69 Vec X, Xold; 70 Vec F, Fold, G; 71 Vec W, Y; 72 73 PetscInt i, k; 74 75 PetscReal fnorm, xnorm = 0, ynorm, gnorm; 76 PetscInt m = qn->m; 77 PetscBool lssucceed; 78 79 PetscScalar rhosc; 80 81 Vec * dX = qn->dX; 82 Vec * dF = qn->dF; 83 PetscScalar * rho = qn->rho; 84 85 /* basically just a regular newton's method except for the application of the jacobian */ 86 PetscFunctionBegin; 87 88 X = snes->vec_sol; /* solution vector */ 89 F = snes->vec_func; /* residual vector */ 90 Y = snes->work[0]; /* work vectors */ 91 G = snes->work[1]; 92 W = snes->work[2]; 93 Xold = snes->work[3]; 94 Fold = snes->work[4]; 95 96 snes->reason = SNES_CONVERGED_ITERATING; 97 98 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 99 snes->iter = 0; 100 snes->norm = 0.; 101 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 102 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 103 if (snes->domainerror) { 104 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 105 PetscFunctionReturn(0); 106 } 107 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 108 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 109 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 110 snes->norm = fnorm; 111 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 112 SNESLogConvHistory(snes,fnorm,0); 113 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 114 115 /* set parameter for default relative tolerance convergence test */ 116 snes->ttol = fnorm*snes->rtol; 117 /* test convergence */ 118 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 119 if (snes->reason) PetscFunctionReturn(0); 120 ierr = VecCopy(F, Fold);CHKERRQ(ierr); 121 ierr = VecCopy(X, Xold);CHKERRQ(ierr); 122 for(i = 0; i < snes->max_its; i++) { 123 /* general purpose update */ 124 if (snes->ops->update) { 125 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 126 } 127 128 /* apply the current iteration of the approximate jacobian */ 129 ierr = LBGFSApplyJinv_Private(snes, i, F, Y);CHKERRQ(ierr); 130 131 /* line search for lambda */ 132 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 133 ynorm = 1; gnorm = fnorm; 134 ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 135 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); 136 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 137 if (snes->domainerror) { 138 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 139 PetscFunctionReturn(0); 140 } 141 if (!lssucceed) { 142 if (++snes->numFailures >= snes->maxFailures) { 143 snes->reason = SNES_DIVERGED_LINE_SEARCH; 144 /* PetscBool ismin; 145 ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 146 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 147 */ 148 break; 149 } 150 } 151 /* Update function and solution vectors */ 152 fnorm = gnorm; 153 ierr = VecCopy(G,F);CHKERRQ(ierr); 154 ierr = VecCopy(W,X);CHKERRQ(ierr); 155 156 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 157 snes->iter = i+1; 158 snes->norm = fnorm; 159 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 160 SNESLogConvHistory(snes,snes->norm,snes->iter); 161 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 162 /* set parameter for default relative tolerance convergence test */ 163 ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 164 if (snes->reason) PetscFunctionReturn(0); 165 166 /* set the differences */ 167 k = i % m; 168 ierr = VecCopy(F, dF[k]);CHKERRQ(ierr); 169 ierr = VecAXPY(dF[k], -1.0, Fold);CHKERRQ(ierr); 170 ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 171 ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 172 ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 173 rho[k] = 1. / rhosc; 174 ierr = VecCopy(F, Fold);CHKERRQ(ierr); 175 ierr = VecCopy(X, Xold);CHKERRQ(ierr); 176 } 177 if (i == snes->max_its) { 178 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 179 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 180 } 181 PetscFunctionReturn(0); 182 } 183 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "SNESSetUp_QN" 187 static PetscErrorCode SNESSetUp_QN(SNES snes) 188 { 189 QNContext * qn = (QNContext *)snes->data; 190 PetscErrorCode ierr; 191 PetscFunctionBegin; 192 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 193 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 194 ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 195 ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr); 196 PetscFunctionReturn(0); 197 } 198 199 #undef __FUNCT__ 200 #define __FUNCT__ "SNESReset_QN" 201 static PetscErrorCode SNESReset_QN(SNES snes) 202 { 203 PetscErrorCode ierr; 204 QNContext * qn; 205 PetscFunctionBegin; 206 if (snes->data) { 207 qn = (QNContext *)snes->data; 208 if (qn->dX) { 209 ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 210 } 211 if (qn->dF) { 212 ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 213 } 214 ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 215 } 216 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 217 PetscFunctionReturn(0); 218 } 219 220 #undef __FUNCT__ 221 #define __FUNCT__ "SNESDestroy_QN" 222 static PetscErrorCode SNESDestroy_QN(SNES snes) 223 { 224 PetscErrorCode ierr; 225 PetscFunctionBegin; 226 ierr = SNESReset_QN(snes);CHKERRQ(ierr); 227 ierr = PetscFree(snes->data);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 #undef __FUNCT__ 232 #define __FUNCT__ "SNESSetFromOptions_QN" 233 static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 234 { 235 236 PetscErrorCode ierr; 237 QNContext * qn; 238 239 PetscFunctionBegin; 240 241 qn = (QNContext *)snes->data; 242 243 ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 244 ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 245 ierr = PetscOptionsTail();CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 249 /* -------------------------------------------------------------------------- */ 250 /*MC 251 SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 252 253 Options Database: 254 255 + -snes_qn_m - Number of past states saved for the L-Broyden methods. 256 257 Notes: This implements the L-BFGS algorithm for the solution of F(x) = 0 using previous change in F(x) and x to 258 form the approximate inverse Jacobian using a series of multiplicative rank-one updates. This will eventually be 259 generalized to implement several limited-memory Broyden methods. 260 261 References: 262 263 L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 264 International Journal of Computer Mathematics, vol. 86, 2009. 265 266 267 Level: beginner 268 269 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 270 271 M*/ 272 EXTERN_C_BEGIN 273 #undef __FUNCT__ 274 #define __FUNCT__ "SNESCreate_QN" 275 PetscErrorCode SNESCreate_QN(SNES snes) 276 { 277 278 PetscErrorCode ierr; 279 QNContext * qn; 280 281 PetscFunctionBegin; 282 snes->ops->setup = SNESSetUp_QN; 283 snes->ops->solve = SNESSolve_QN; 284 snes->ops->destroy = SNESDestroy_QN; 285 snes->ops->setfromoptions = SNESSetFromOptions_QN; 286 snes->ops->view = 0; 287 snes->ops->reset = SNESReset_QN; 288 289 snes->usespc = PETSC_TRUE; 290 snes->usesksp = PETSC_FALSE; 291 292 ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 293 snes->data = (void *) qn; 294 qn->m = 100; 295 qn->dX = PETSC_NULL; 296 qn->dF = PETSC_NULL; 297 298 snes->ops->linesearchcubic = SNESLineSearchCubic; 299 snes->ops->linesearchquadratic = SNESLineSearchQuadratic; 300 snes->ops->linesearchno = SNESLineSearchNo; 301 snes->ops->linesearchnonorms = SNESLineSearchNoNorms; 302 303 ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr); 304 305 PetscFunctionReturn(0); 306 } 307 EXTERN_C_END 308