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