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