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