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