1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 3 #define H(i,j) qn->dXdFmat[i*qn->m + j] 4 5 const char *const SNESQNScaleTypes[] = {"NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0}; 6 const char *const SNESQNRestartTypes[] = {"NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0}; 7 const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_RESTART_",0}; 8 9 typedef enum {SNES_QN_LBFGS = 0, 10 SNES_QN_BROYDEN = 1, 11 SNES_QN_BADBROYDEN = 2 12 } SNESQNType; 13 14 typedef struct { 15 Vec *U; /* Stored past states (vary from method to method) */ 16 Vec *V; /* Stored past states (vary from method to method) */ 17 PetscInt m; /* The number of kept previous steps */ 18 PetscScalar *alpha, *beta; 19 PetscScalar *dXtdF, *dFtdX, *YtdX; 20 PetscBool singlereduction; /* Aggregated reduction implementation */ 21 PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */ 22 PetscViewer monitor; 23 PetscReal powell_gamma; /* Powell angle restart condition */ 24 PetscReal powell_downhill; /* Powell descent restart condition */ 25 PetscReal scaling; /* scaling of H0 */ 26 SNESQNType type; /* the type of quasi-newton method used */ 27 SNESQNScaleType scale_type; /* the type of scaling used */ 28 SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */ 29 } SNES_QN; 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "SNESQNApply_Broyden" 33 PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold, Vec D, Vec Dold) 34 { 35 PetscErrorCode ierr; 36 SNES_QN *qn = (SNES_QN*)snes->data; 37 Vec W = snes->work[3]; 38 Vec *U = qn->U; 39 Vec *V = qn->V; 40 KSPConvergedReason kspreason; 41 PetscInt k,i,lits; 42 PetscInt m = qn->m; 43 PetscScalar gdot; 44 PetscInt l = m; 45 46 PetscFunctionBegin; 47 if (it < m) l = it; 48 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 49 ierr = KSPSolve(snes->ksp,D,W);CHKERRQ(ierr); 50 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 51 if (kspreason < 0) { 52 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 53 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 54 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 55 PetscFunctionReturn(0); 56 } 57 } 58 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 59 snes->linear_its += lits; 60 ierr = VecCopy(W,Y);CHKERRQ(ierr); 61 } else { 62 ierr = VecCopy(D,Y);CHKERRQ(ierr); 63 ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 64 } 65 66 /* inward recursion starting at the first update and working forward */ 67 if (it > 1) { 68 for (i = 0; i < l-1; i++) { 69 k = (it+i-l)%l; 70 ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 71 ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 72 if (qn->monitor) { 73 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 74 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 75 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 76 } 77 } 78 } 79 if (it < m) l = it; 80 81 /* set W to be the last step, post-linesearch */ 82 ierr = VecCopy(Xold,W);CHKERRQ(ierr); 83 ierr = VecAXPY(W,-1.0,X);CHKERRQ(ierr); 84 85 if (l > 0) { 86 k = (it-1)%l; 87 ierr = VecCopy(W,U[k]);CHKERRQ(ierr); 88 ierr = VecAXPY(W,-1.0,Y);CHKERRQ(ierr); 89 ierr = VecDot(U[k],W,&gdot);CHKERRQ(ierr); 90 ierr = VecCopy(Y,V[k]);CHKERRQ(ierr); 91 ierr = VecScale(V[k],1.0/gdot);CHKERRQ(ierr); 92 ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 93 ierr = VecAXPY(Y,gdot,V[k]);CHKERRQ(ierr); 94 if (qn->monitor) { 95 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 96 ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 97 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 98 } 99 } 100 PetscFunctionReturn(0); 101 } 102 103 #undef __FUNCT__ 104 #define __FUNCT__ "SNESQNApply_BadBroyden" 105 PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 106 { 107 PetscErrorCode ierr; 108 SNES_QN *qn = (SNES_QN*)snes->data; 109 Vec W = snes->work[3]; 110 Vec *U = qn->U; 111 Vec *T = qn->V; 112 113 /* ksp thing for jacobian scaling */ 114 KSPConvergedReason kspreason; 115 PetscInt k, i, lits; 116 PetscInt m = qn->m; 117 PetscScalar gdot; 118 PetscInt l = m; 119 120 PetscFunctionBegin; 121 if (it < m) l = it; 122 ierr = VecCopy(D,Y);CHKERRQ(ierr); 123 if (l > 0) { 124 k = (it-1)%l; 125 ierr = VecCopy(Dold,U[k]);CHKERRQ(ierr); 126 ierr = VecAXPY(U[k],-1.0,D);CHKERRQ(ierr); 127 ierr = VecDot(U[k],U[k],&gdot);CHKERRQ(ierr); 128 ierr = VecCopy(D,T[k]);CHKERRQ(ierr); 129 ierr = VecScale(T[k],1.0/gdot);CHKERRQ(ierr); 130 ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 131 ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 132 if (qn->monitor) { 133 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 134 ierr = PetscViewerASCIIPrintf(qn->monitor, "update: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 135 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 136 } 137 } 138 139 /* inward recursion starting at the first update and working forward */ 140 if (it > 1) { 141 for (i = 0; i < l-1; i++) { 142 k = (it+i-l)%l; 143 ierr = VecDot(U[k],Y,&gdot);CHKERRQ(ierr); 144 ierr = VecAXPY(Y,gdot,T[k]);CHKERRQ(ierr); 145 if (qn->monitor) { 146 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 147 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));CHKERRQ(ierr); 148 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 149 } 150 } 151 } 152 153 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 154 ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 155 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 156 if (kspreason < 0) { 157 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 158 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 159 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 160 PetscFunctionReturn(0); 161 } 162 } 163 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 164 snes->linear_its += lits; 165 ierr = VecCopy(W,Y);CHKERRQ(ierr); 166 } else { 167 ierr = VecScale(Y,qn->scaling);CHKERRQ(ierr); 168 } 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #define __FUNCT__ "SNESQNApply_LBFGS" 174 PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) 175 { 176 PetscErrorCode ierr; 177 SNES_QN *qn = (SNES_QN*)snes->data; 178 Vec W = snes->work[3]; 179 Vec *dX = qn->U; 180 Vec *dF = qn->V; 181 PetscScalar *alpha = qn->alpha; 182 PetscScalar *beta = qn->beta; 183 PetscScalar *dXtdF = qn->dXtdF; 184 PetscScalar *dFtdX = qn->dFtdX; 185 PetscScalar *YtdX = qn->YtdX; 186 187 /* ksp thing for jacobian scaling */ 188 KSPConvergedReason kspreason; 189 PetscInt k,i,j,g,lits; 190 PetscInt m = qn->m; 191 PetscScalar t; 192 PetscInt l = m; 193 194 PetscFunctionBegin; 195 if (it < m) l = it; 196 if (it > 0) { 197 k = (it - 1) % l; 198 ierr = VecCopy(D, dF[k]);CHKERRQ(ierr); 199 ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); 200 ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); 201 ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); 202 if (qn->singlereduction) { 203 ierr = VecMDot(dF[k], l, dX, dXtdF);CHKERRQ(ierr); 204 ierr = VecMDot(dX[k], l, dF, dFtdX);CHKERRQ(ierr); 205 for (j = 0; j < l; j++) { 206 H(k, j) = dFtdX[j]; 207 H(j, k) = dXtdF[j]; 208 } 209 /* copy back over to make the computation of alpha and beta easier */ 210 for (j = 0; j < l; j++) dXtdF[j] = H(j, j); 211 } else { 212 ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); 213 } 214 if (qn->scale_type == SNES_QN_SCALE_SHANNO) { 215 PetscReal dFtdF; 216 ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); 217 qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; 218 } else if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 219 ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); 220 } 221 } 222 223 ierr = VecCopy(D,Y);CHKERRQ(ierr); 224 225 if (qn->singlereduction) { 226 ierr = VecMDot(Y,l,dX,YtdX);CHKERRQ(ierr); 227 } 228 /* outward recursion starting at iteration k's update and working back */ 229 for (i=0; i<l; i++) { 230 k = (it-i-1)%l; 231 if (qn->singlereduction) { 232 /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ 233 t = YtdX[k]; 234 for (j=0; j<i; j++) { 235 g = (it-j-1)%l; 236 t += -alpha[g]*H(g, k); 237 } 238 alpha[k] = t/H(k,k); 239 } else { 240 ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); 241 alpha[k] = t/dXtdF[k]; 242 } 243 if (qn->monitor) { 244 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 245 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); 246 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 247 } 248 ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); 249 } 250 251 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 252 ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); 253 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 254 if (kspreason < 0) { 255 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 256 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 257 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 258 PetscFunctionReturn(0); 259 } 260 } 261 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 262 snes->linear_its += lits; 263 ierr = VecCopy(W, Y);CHKERRQ(ierr); 264 } else { 265 ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); 266 } 267 if (qn->singlereduction) { 268 ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); 269 } 270 /* inward recursion starting at the first update and working forward */ 271 for (i = 0; i < l; i++) { 272 k = (it + i - l) % l; 273 if (qn->singlereduction) { 274 t = YtdX[k]; 275 for (j = 0; j < i; j++) { 276 g = (it + j - l) % l; 277 t += (alpha[g] - beta[g])*H(k, g); 278 } 279 beta[k] = t / H(k, k); 280 } else { 281 ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); 282 beta[k] = t / dXtdF[k]; 283 } 284 ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); 285 if (qn->monitor) { 286 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 287 ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); 288 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 289 } 290 } 291 PetscFunctionReturn(0); 292 } 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "SNESSolve_QN" 296 static PetscErrorCode SNESSolve_QN(SNES snes) 297 { 298 PetscErrorCode ierr; 299 SNES_QN *qn = (SNES_QN*) snes->data; 300 Vec X,Xold; 301 Vec F; 302 Vec Y,D,Dold; 303 PetscInt i, i_r; 304 PetscReal fnorm,xnorm,ynorm,gnorm; 305 PetscBool lssucceed,powell,periodic; 306 PetscScalar DolddotD,DolddotDold,DdotD,YdotD; 307 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 308 SNESConvergedReason reason; 309 310 /* basically just a regular newton's method except for the application of the jacobian */ 311 312 PetscFunctionBegin; 313 F = snes->vec_func; /* residual vector */ 314 Y = snes->vec_sol_update; /* search direction generated by J^-1D*/ 315 316 X = snes->vec_sol; /* solution vector */ 317 Xold = snes->work[0]; 318 319 /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */ 320 D = snes->work[1]; 321 Dold = snes->work[2]; 322 323 snes->reason = SNES_CONVERGED_ITERATING; 324 325 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 326 snes->iter = 0; 327 snes->norm = 0.; 328 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 329 330 if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 331 ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 332 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 333 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 334 snes->reason = SNES_DIVERGED_INNER; 335 PetscFunctionReturn(0); 336 } 337 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 338 } else { 339 if (!snes->vec_func_init_set) { 340 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 341 if (snes->domainerror) { 342 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 343 PetscFunctionReturn(0); 344 } 345 } else snes->vec_func_init_set = PETSC_FALSE; 346 if (!snes->norm_init_set) { 347 /* convergence test */ 348 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 349 if (PetscIsInfOrNanReal(fnorm)) { 350 snes->reason = SNES_DIVERGED_FNORM_NAN; 351 PetscFunctionReturn(0); 352 } 353 } else { 354 fnorm = snes->norm_init; 355 snes->norm_init_set = PETSC_FALSE; 356 } 357 } 358 if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 359 ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 360 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 361 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 362 snes->reason = SNES_DIVERGED_INNER; 363 PetscFunctionReturn(0); 364 } 365 } else { 366 ierr = VecCopy(F,D);CHKERRQ(ierr); 367 } 368 ierr = VecCopy(D,Y);CHKERRQ(ierr); 369 370 371 ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 372 snes->norm = fnorm; 373 ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 374 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 375 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 376 377 /* set parameter for default relative tolerance convergence test */ 378 snes->ttol = fnorm*snes->rtol; 379 /* test convergence */ 380 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 381 if (snes->reason) PetscFunctionReturn(0); 382 383 /* scale the initial update */ 384 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 385 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 386 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 387 } 388 389 for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) { 390 switch (qn->type) { 391 case SNES_QN_BADBROYDEN: 392 ierr = SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 393 break; 394 case SNES_QN_BROYDEN: 395 ierr = SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 396 break; 397 case SNES_QN_LBFGS: 398 SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);CHKERRQ(ierr); 399 break; 400 } 401 /* line search for lambda */ 402 ynorm = 1; gnorm = fnorm; 403 ierr = VecCopy(D, Dold);CHKERRQ(ierr); 404 ierr = VecCopy(X, Xold);CHKERRQ(ierr); 405 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 406 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 407 if (snes->domainerror) { 408 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 409 PetscFunctionReturn(0); 410 } 411 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssucceed);CHKERRQ(ierr); 412 if (!lssucceed) { 413 if (++snes->numFailures >= snes->maxFailures) { 414 snes->reason = SNES_DIVERGED_LINE_SEARCH; 415 break; 416 } 417 } 418 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 419 if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { 420 ierr = SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);CHKERRQ(ierr); 421 } 422 423 /* convergence monitoring */ 424 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); 425 426 ierr = SNESSetIterationNumber(snes, i+1);CHKERRQ(ierr); 427 ierr = SNESSetFunctionNorm(snes, fnorm);CHKERRQ(ierr); 428 429 ierr = SNESLogConvergenceHistory(snes,snes->norm,snes->iter);CHKERRQ(ierr); 430 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 431 /* set parameter for default relative tolerance convergence test */ 432 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 433 if (snes->reason) PetscFunctionReturn(0); 434 if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 435 ierr = SNESApplyPC(snes,X,F,&fnorm,D);CHKERRQ(ierr); 436 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 437 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 438 snes->reason = SNES_DIVERGED_INNER; 439 PetscFunctionReturn(0); 440 } 441 } else { 442 ierr = VecCopy(F, D);CHKERRQ(ierr); 443 } 444 445 powell = PETSC_FALSE; 446 if (qn->restart_type == SNES_QN_RESTART_POWELL) { 447 /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */ 448 ierr = VecDotBegin(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 449 ierr = VecDotBegin(Dold, D, &DolddotD);CHKERRQ(ierr); 450 ierr = VecDotBegin(D, D, &DdotD);CHKERRQ(ierr); 451 ierr = VecDotBegin(Y, D, &YdotD);CHKERRQ(ierr); 452 ierr = VecDotEnd(Dold, Dold, &DolddotDold);CHKERRQ(ierr); 453 ierr = VecDotEnd(Dold, D, &DolddotD);CHKERRQ(ierr); 454 ierr = VecDotEnd(D, D, &DdotD);CHKERRQ(ierr); 455 ierr = VecDotEnd(Y, D, &YdotD);CHKERRQ(ierr); 456 if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE; 457 } 458 periodic = PETSC_FALSE; 459 if (qn->restart_type == SNES_QN_RESTART_PERIODIC) { 460 if (i_r>qn->m-1) periodic = PETSC_TRUE; 461 } 462 /* restart if either powell or periodic restart is satisfied. */ 463 if (powell || periodic) { 464 if (qn->monitor) { 465 ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 466 ierr = PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);CHKERRQ(ierr); 467 ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); 468 } 469 i_r = -1; 470 /* general purpose update */ 471 if (snes->ops->update) { 472 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 473 } 474 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 475 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 476 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 477 } 478 } 479 /* general purpose update */ 480 if (snes->ops->update) { 481 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 482 } 483 } 484 if (i == snes->max_its) { 485 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 486 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 487 } 488 PetscFunctionReturn(0); 489 } 490 491 #undef __FUNCT__ 492 #define __FUNCT__ "SNESSetUp_QN" 493 static PetscErrorCode SNESSetUp_QN(SNES snes) 494 { 495 SNES_QN *qn = (SNES_QN*)snes->data; 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESQN only supports left preconditioning"); 500 501 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);CHKERRQ(ierr); 502 ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);CHKERRQ(ierr); 503 ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, 504 qn->m, PetscScalar, &qn->beta, 505 qn->m, PetscScalar, &qn->dXtdF);CHKERRQ(ierr); 506 507 if (qn->singlereduction) { 508 ierr = PetscMalloc3(qn->m*qn->m, PetscScalar, &qn->dXdFmat, 509 qn->m, PetscScalar, &qn->dFtdX, 510 qn->m, PetscScalar, &qn->YtdX);CHKERRQ(ierr); 511 } 512 ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 513 514 /* set up the line search */ 515 if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { 516 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 517 } 518 PetscFunctionReturn(0); 519 } 520 521 #undef __FUNCT__ 522 #define __FUNCT__ "SNESReset_QN" 523 static PetscErrorCode SNESReset_QN(SNES snes) 524 { 525 PetscErrorCode ierr; 526 SNES_QN *qn; 527 528 PetscFunctionBegin; 529 if (snes->data) { 530 qn = (SNES_QN*)snes->data; 531 if (qn->U) { 532 ierr = VecDestroyVecs(qn->m, &qn->U);CHKERRQ(ierr); 533 } 534 if (qn->V) { 535 ierr = VecDestroyVecs(qn->m, &qn->V);CHKERRQ(ierr); 536 } 537 if (qn->singlereduction) { 538 ierr = PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);CHKERRQ(ierr); 539 } 540 ierr = PetscFree3(qn->alpha, qn->beta, qn->dXtdF);CHKERRQ(ierr); 541 } 542 PetscFunctionReturn(0); 543 } 544 545 #undef __FUNCT__ 546 #define __FUNCT__ "SNESDestroy_QN" 547 static PetscErrorCode SNESDestroy_QN(SNES snes) 548 { 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 ierr = SNESReset_QN(snes);CHKERRQ(ierr); 553 ierr = PetscFree(snes->data);CHKERRQ(ierr); 554 ierr = PetscObjectComposeFunction((PetscObject)snes,"",NULL);CHKERRQ(ierr); 555 PetscFunctionReturn(0); 556 } 557 558 #undef __FUNCT__ 559 #define __FUNCT__ "SNESSetFromOptions_QN" 560 static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 561 { 562 563 PetscErrorCode ierr; 564 SNES_QN *qn = (SNES_QN*)snes->data; 565 PetscBool monflg = PETSC_FALSE,flg; 566 SNESLineSearch linesearch; 567 SNESQNRestartType rtype = qn->restart_type; 568 SNESQNScaleType stype = qn->scale_type; 569 570 PetscFunctionBegin; 571 ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 572 ierr = PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);CHKERRQ(ierr); 573 ierr = PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);CHKERRQ(ierr); 574 ierr = PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);CHKERRQ(ierr); 575 ierr = PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);CHKERRQ(ierr); 576 ierr = PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);CHKERRQ(ierr); 577 ierr = PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);CHKERRQ(ierr); 578 if (flg) ierr = SNESQNSetScaleType(snes,stype);CHKERRQ(ierr); 579 580 ierr = PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);CHKERRQ(ierr); 581 if (flg) ierr = SNESQNSetRestartType(snes,rtype);CHKERRQ(ierr); 582 583 ierr = PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes, 584 (PetscEnum)qn->type,(PetscEnum*)&qn->type,NULL);CHKERRQ(ierr); 585 ierr = PetscOptionsTail();CHKERRQ(ierr); 586 if (!snes->linesearch) { 587 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 588 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 589 } 590 if (monflg) { 591 qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 592 } 593 PetscFunctionReturn(0); 594 } 595 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "SNESQNSetRestartType" 599 /*@ 600 SNESQNSetRestartType - Sets the restart type for SNESQN. 601 602 Logically Collective on SNES 603 604 Input Parameters: 605 + snes - the iterative context 606 - rtype - restart type 607 608 Options Database: 609 + -snes_qn_restart_type<powell,periodic,none> - set the restart type 610 - -snes_qn_restart[10] - sets the number of iterations before restart for periodic 611 612 Level: intermediate 613 614 SNESQNRestartTypes: 615 + SNES_QN_RESTART_NONE - never restart 616 . SNES_QN_RESTART_POWELL - restart based upon descent criteria 617 - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations 618 619 Notes: 620 The default line search used is the L2 line search and it requires two additional function evaluations. 621 622 .keywords: SNES, SNESQN, restart, type, set SNESLineSearch 623 @*/ 624 PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype) 625 { 626 PetscErrorCode ierr; 627 628 PetscFunctionBegin; 629 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 630 ierr = PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));CHKERRQ(ierr); 631 PetscFunctionReturn(0); 632 } 633 634 #undef __FUNCT__ 635 #define __FUNCT__ "SNESQNSetScaleType" 636 /*@ 637 SNESQNSetScaleType - Sets the scaling type for the inner inverse jacobian in SNESQN. 638 639 Logically Collective on SNES 640 641 Input Parameters: 642 + snes - the iterative context 643 - stype - scale type 644 645 Options Database: 646 . -snes_qn_scale_type<shanno,none,linesearch,jacobian> 647 648 Level: intermediate 649 650 SNESQNSelectTypes: 651 + SNES_QN_SCALE_NONE - don't scale the problem 652 . SNES_QN_SCALE_SHANNO - use shanno scaling 653 . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda 654 - SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian. 655 656 .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch 657 @*/ 658 659 PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype) 660 { 661 PetscErrorCode ierr; 662 663 PetscFunctionBegin; 664 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 665 ierr = PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));CHKERRQ(ierr); 666 PetscFunctionReturn(0); 667 } 668 669 #undef __FUNCT__ 670 #define __FUNCT__ "SNESQNSetScaleType_QN" 671 PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype) 672 { 673 SNES_QN *qn = (SNES_QN*)snes->data; 674 675 PetscFunctionBegin; 676 qn->scale_type = stype; 677 PetscFunctionReturn(0); 678 } 679 680 #undef __FUNCT__ 681 #define __FUNCT__ "SNESQNSetRestartType_QN" 682 PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype) 683 { 684 SNES_QN *qn = (SNES_QN*)snes->data; 685 686 PetscFunctionBegin; 687 qn->restart_type = rtype; 688 PetscFunctionReturn(0); 689 } 690 691 /* -------------------------------------------------------------------------- */ 692 /*MC 693 SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 694 695 Options Database: 696 697 + -snes_qn_m - Number of past states saved for the L-Broyden methods. 698 . -snes_qn_powell_angle - Angle condition for restart. 699 . -snes_qn_powell_descent - Descent condition for restart. 700 . -snes_linesearch_type <cp, l2, basic> - Type of line search. 701 - -snes_qn_monitor - Monitors the quasi-newton jacobian. 702 703 Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using 704 previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one 705 updates. 706 707 When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of 708 these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this 709 iteration as the current iteration's values when constructing the approximate jacobian. The second, composed, 710 perturbs the problem the jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner. 711 712 References: 713 714 L-Broyden Methods: a generalization of the L-BFGS method to the limited memory Broyden family, M. B. Reed, 715 International Journal of Computer Mathematics, vol. 86, 2009. 716 717 Experiments with Quasi-Newton Methods in Solving Stiff ODE Systems, Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker 718 SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985. 719 720 721 Level: beginner 722 723 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR 724 725 M*/ 726 #undef __FUNCT__ 727 #define __FUNCT__ "SNESCreate_QN" 728 PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes) 729 { 730 PetscErrorCode ierr; 731 SNES_QN *qn; 732 733 PetscFunctionBegin; 734 snes->ops->setup = SNESSetUp_QN; 735 snes->ops->solve = SNESSolve_QN; 736 snes->ops->destroy = SNESDestroy_QN; 737 snes->ops->setfromoptions = SNESSetFromOptions_QN; 738 snes->ops->view = 0; 739 snes->ops->reset = SNESReset_QN; 740 741 snes->pcside = PC_LEFT; 742 snes->functype = SNES_FUNCTION_PRECONDITIONED; 743 744 snes->usespc = PETSC_TRUE; 745 snes->usesksp = PETSC_FALSE; 746 747 if (!snes->tolerancesset) { 748 snes->max_funcs = 30000; 749 snes->max_its = 10000; 750 } 751 752 ierr = PetscNewLog(snes,SNES_QN,&qn);CHKERRQ(ierr); 753 snes->data = (void*) qn; 754 qn->m = 10; 755 qn->scaling = 1.0; 756 qn->U = NULL; 757 qn->V = NULL; 758 qn->dXtdF = NULL; 759 qn->dFtdX = NULL; 760 qn->dXdFmat = NULL; 761 qn->monitor = NULL; 762 qn->singlereduction = PETSC_FALSE; 763 qn->powell_gamma = 0.9999; 764 qn->scale_type = SNES_QN_SCALE_SHANNO; 765 qn->restart_type = SNES_QN_RESTART_POWELL; 766 qn->type = SNES_QN_LBFGS; 767 768 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);CHKERRQ(ierr); 769 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);CHKERRQ(ierr); 770 PetscFunctionReturn(0); 771 } 772