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