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