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