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