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