1 #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/ 2 const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0}; 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESReset_NCG" 6 PetscErrorCode SNESReset_NCG(SNES snes) 7 { 8 9 PetscFunctionBegin; 10 PetscFunctionReturn(0); 11 } 12 13 #define SNESLINESEARCHNCGLINEAR "linear" 14 15 /* 16 SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 17 18 Input Parameter: 19 . snes - the SNES context 20 21 Application Interface Routine: SNESDestroy() 22 */ 23 #undef __FUNCT__ 24 #define __FUNCT__ "SNESDestroy_NCG" 25 PetscErrorCode SNESDestroy_NCG(SNES snes) 26 { 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 ierr = PetscFree(snes->data);CHKERRQ(ierr); 31 PetscFunctionReturn(0); 32 } 33 34 /* 35 SNESSetUp_NCG - Sets up the internal data structures for the later use 36 of the SNESNCG nonlinear solver. 37 38 Input Parameters: 39 + snes - the SNES context 40 - x - the solution vector 41 42 Application Interface Routine: SNESSetUp() 43 */ 44 45 EXTERN_C_BEGIN 46 extern PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch); 47 EXTERN_C_END 48 49 #undef __FUNCT__ 50 #define __FUNCT__ "SNESSetUp_NCG" 51 PetscErrorCode SNESSetUp_NCG(SNES snes) 52 { 53 PetscErrorCode ierr; 54 55 PetscFunctionBegin; 56 ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr); 57 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 58 ierr = SNESLineSearchRegisterDynamic(SNESLINESEARCHNCGLINEAR, PETSC_NULL,"SNESLineSearchCreate_NCGLinear", SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr); 59 60 PetscFunctionReturn(0); 61 } 62 /* 63 SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 64 65 Input Parameter: 66 . snes - the SNES context 67 68 Application Interface Routine: SNESSetFromOptions() 69 */ 70 #undef __FUNCT__ 71 #define __FUNCT__ "SNESSetFromOptions_NCG" 72 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 73 { 74 SNES_NCG *ncg = (SNES_NCG *)snes->data; 75 PetscErrorCode ierr; 76 PetscBool debug; 77 SNESLineSearch linesearch; 78 SNESNCGType ncgtype=ncg->type; 79 80 PetscFunctionBegin; 81 ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 82 ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 83 ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,PETSC_NULL);CHKERRQ(ierr); 84 ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 85 if (debug) { 86 ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 87 } 88 ierr = PetscOptionsTail();CHKERRQ(ierr); 89 if (!snes->linesearch) { 90 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 91 if (!snes->pc) { 92 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 93 } else { 94 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 95 } 96 } 97 PetscFunctionReturn(0); 98 } 99 100 /* 101 SNESView_NCG - Prints info from the SNESNCG data structure. 102 103 Input Parameters: 104 + SNES - the SNES context 105 - viewer - visualization context 106 107 Application Interface Routine: SNESView() 108 */ 109 #undef __FUNCT__ 110 #define __FUNCT__ "SNESView_NCG" 111 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 112 { 113 PetscBool iascii; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 118 if (iascii) { 119 } 120 PetscFunctionReturn(0); 121 } 122 123 #undef __FUNCT__ 124 #define __FUNCT__ "SNESLineSearchApply_NCGLinear" 125 PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 126 { 127 PetscScalar alpha, ptAp; 128 Vec X, Y, F, W; 129 SNES snes; 130 PetscErrorCode ierr; 131 PetscReal *fnorm, *xnorm, *ynorm; 132 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 133 134 PetscFunctionBegin; 135 136 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 137 X = linesearch->vec_sol; 138 W = linesearch->vec_sol_new; 139 F = linesearch->vec_func; 140 Y = linesearch->vec_update; 141 fnorm = &linesearch->fnorm; 142 xnorm = &linesearch->xnorm; 143 ynorm = &linesearch->ynorm; 144 145 /* 146 147 The exact step size for unpreconditioned linear CG is just: 148 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 149 */ 150 ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr); 151 ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 152 ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 153 ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 154 alpha = alpha / ptAp; 155 ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr); 156 ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr); 157 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 158 159 ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 160 ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr); 161 ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr); 162 PetscFunctionReturn(0); 163 } 164 165 EXTERN_C_BEGIN 166 #undef __FUNCT__ 167 #define __FUNCT__ "SNESLineSearchCreate_NCGLinear" 168 169 PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 170 { 171 PetscFunctionBegin; 172 linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 173 linesearch->ops->destroy = PETSC_NULL; 174 linesearch->ops->setfromoptions = PETSC_NULL; 175 linesearch->ops->reset = PETSC_NULL; 176 linesearch->ops->view = PETSC_NULL; 177 linesearch->ops->setup = PETSC_NULL; 178 PetscFunctionReturn(0); 179 } 180 EXTERN_C_END 181 182 #undef __FUNCT__ 183 #define __FUNCT__ "SNESNCGComputeYtJtF_Private" 184 /* 185 186 Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 187 188 */ 189 PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) 190 { 191 PetscErrorCode ierr; 192 PetscScalar ftf, ftg, fty, h; 193 194 PetscFunctionBegin; 195 ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 196 ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 197 h = 1e-5*fty / fty; 198 ierr = VecCopy(X, W);CHKERRQ(ierr); 199 ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 200 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 201 ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 202 *ytJtf = (ftg - ftf) / h; 203 PetscFunctionReturn(0); 204 } 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "SNESNCGSetType" 208 /*@ 209 SNESNCGSetType - Sets the conjugate update type for SNESNCG. 210 211 Logically Collective on SNES 212 213 Input Parameters: 214 + snes - the iterative context 215 - btype - update type 216 217 Options Database: 218 . -snes_ncg_type<prp,fr,hs,dy,cd> 219 220 Level: intermediate 221 222 SNESNCGSelectTypes: 223 + SNES_NCG_FR - Fletcher-Reeves update 224 . SNES_NCG_PRP - Polak-Ribiere-Polyak update 225 . SNES_NCG_HS - Hestenes-Steifel update 226 . SNES_NCG_DY - Dai-Yuan update 227 - SNES_NCG_CD - Conjugate Descent update 228 229 Notes: 230 PRP is the default, and the only one that tolerates generalized search directions. 231 232 .keywords: SNES, SNESNCG, selection, type, set 233 @*/ 234 PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) { 235 PetscErrorCode ierr; 236 PetscFunctionBegin; 237 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 238 ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr); 239 PetscFunctionReturn(0); 240 } 241 242 243 EXTERN_C_BEGIN 244 #undef __FUNCT__ 245 #define __FUNCT__ "SNESNCGSetType_NCG" 246 PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) { 247 SNES_NCG *ncg = (SNES_NCG *)snes->data; 248 PetscFunctionBegin; 249 ncg->type = btype; 250 PetscFunctionReturn(0); 251 } 252 EXTERN_C_END 253 254 /* 255 SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 256 257 Input Parameters: 258 . snes - the SNES context 259 260 Output Parameter: 261 . outits - number of iterations until termination 262 263 Application Interface Routine: SNESSolve() 264 */ 265 #undef __FUNCT__ 266 #define __FUNCT__ "SNESSolve_NCG" 267 PetscErrorCode SNESSolve_NCG(SNES snes) 268 { 269 SNES_NCG *ncg = (SNES_NCG *)snes->data; 270 Vec X, dX, lX, F, B, Fold; 271 PetscReal fnorm, ynorm, xnorm, beta = 0.0; 272 PetscScalar dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold; 273 PetscInt maxits, i; 274 PetscErrorCode ierr; 275 SNESConvergedReason reason; 276 PetscBool lsSuccess = PETSC_TRUE; 277 SNESLineSearch linesearch; 278 279 PetscFunctionBegin; 280 snes->reason = SNES_CONVERGED_ITERATING; 281 282 maxits = snes->max_its; /* maximum number of iterations */ 283 X = snes->vec_sol; /* X^n */ 284 Fold = snes->work[0]; /* The previous iterate of X */ 285 dX = snes->work[1]; /* the preconditioned direction */ 286 lX = snes->vec_sol_update; /* the conjugate direction */ 287 F = snes->vec_func; /* residual vector */ 288 B = snes->vec_rhs; /* the right hand side */ 289 290 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 291 292 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 293 snes->iter = 0; 294 snes->norm = 0.; 295 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 296 297 /* compute the initial function and preconditioned update dX */ 298 if (!snes->vec_func_init_set) { 299 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 300 if (snes->domainerror) { 301 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 302 PetscFunctionReturn(0); 303 } 304 } else { 305 snes->vec_func_init_set = PETSC_FALSE; 306 } 307 if (!snes->norm_init_set) { 308 /* convergence test */ 309 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 310 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 311 } else { 312 fnorm = snes->norm_init; 313 snes->norm_init_set = PETSC_FALSE; 314 } 315 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 316 snes->norm = fnorm; 317 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 318 SNESLogConvHistory(snes,fnorm,0); 319 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 320 321 /* set parameter for default relative tolerance convergence test */ 322 snes->ttol = fnorm*snes->rtol; 323 /* test convergence */ 324 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 325 if (snes->reason) PetscFunctionReturn(0); 326 327 /* Call general purpose update function */ 328 if (snes->ops->update) { 329 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 330 } 331 332 /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 333 334 if (snes->pc && snes->pcside == PC_RIGHT) { 335 ierr = VecCopy(X, dX);CHKERRQ(ierr); 336 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 337 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 338 ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr); 339 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 340 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 341 snes->reason = SNES_DIVERGED_INNER; 342 PetscFunctionReturn(0); 343 } 344 ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr); 345 } else { 346 ierr = VecCopy(F, dX);CHKERRQ(ierr); 347 } 348 ierr = VecCopy(dX, lX);CHKERRQ(ierr); 349 ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr); 350 /* 351 } else { 352 ierr = SNESNCGComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr); 353 } 354 */ 355 for (i = 1; i < maxits + 1; i++) { 356 lsSuccess = PETSC_TRUE; 357 /* some update types require the old update direction or conjugate direction */ 358 if (ncg->type != SNES_NCG_FR) { 359 ierr = VecCopy(F, Fold);CHKERRQ(ierr); 360 } 361 ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, lX);CHKERRQ(ierr); 362 ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr); 363 if (!lsSuccess) { 364 if (++snes->numFailures >= snes->maxFailures) { 365 snes->reason = SNES_DIVERGED_LINE_SEARCH; 366 PetscFunctionReturn(0); 367 } 368 } 369 if (snes->nfuncs >= snes->max_funcs) { 370 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 371 PetscFunctionReturn(0); 372 } 373 if (snes->domainerror) { 374 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 375 PetscFunctionReturn(0); 376 } 377 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 378 /* Monitor convergence */ 379 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 380 snes->iter = i; 381 snes->norm = fnorm; 382 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 383 SNESLogConvHistory(snes,snes->norm,0); 384 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 385 386 /* Test for convergence */ 387 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 388 if (snes->reason) PetscFunctionReturn(0); 389 390 /* Call general purpose update function */ 391 if (snes->ops->update) { 392 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 393 } 394 if (snes->pc && snes->pcside == PC_RIGHT) { 395 ierr = VecCopy(X,dX);CHKERRQ(ierr); 396 ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 397 ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 398 ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr); 399 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 400 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 401 snes->reason = SNES_DIVERGED_INNER; 402 PetscFunctionReturn(0); 403 } 404 ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr); 405 } else { 406 ierr = VecCopy(F, dX);CHKERRQ(ierr); 407 } 408 409 /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 410 switch(ncg->type) { 411 case SNES_NCG_FR: /* Fletcher-Reeves */ 412 dXolddotFold = dXdotF; 413 ierr = VecDot(dX, dX, &dXdotF);CHKERRQ(ierr); 414 beta = PetscRealPart(dXdotF / dXolddotFold); 415 break; 416 case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 417 dXolddotFold = dXdotF; 418 ierr = VecDotBegin(F, dX, &dXdotF);CHKERRQ(ierr); 419 ierr = VecDotBegin(Fold, dX, &dXdotFold);CHKERRQ(ierr); 420 ierr = VecDotEnd(F, dX, &dXdotF);CHKERRQ(ierr); 421 ierr = VecDotEnd(Fold, dX, &dXdotFold);CHKERRQ(ierr); 422 beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold)); 423 if (beta < 0.0) beta = 0.0; /* restart */ 424 break; 425 case SNES_NCG_HS: /* Hestenes-Stiefel */ 426 ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr); 427 ierr = VecDotBegin(dX, Fold, &dXdotFold);CHKERRQ(ierr); 428 ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr); 429 ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr); 430 ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr); 431 ierr = VecDotEnd(dX, Fold, &dXdotFold);CHKERRQ(ierr); 432 ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr); 433 ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr); 434 beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold)); 435 break; 436 case SNES_NCG_DY: /* Dai-Yuan */ 437 ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr); 438 ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr); 439 ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr); 440 ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr); 441 ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr); 442 ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr); 443 beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr); 444 break; 445 case SNES_NCG_CD: /* Conjugate Descent */ 446 ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr); 447 ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr); 448 ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr); 449 ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr); 450 beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr); 451 break; 452 } 453 if (ncg->monitor) { 454 ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr); 455 } 456 ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 457 } 458 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 459 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 460 PetscFunctionReturn(0); 461 } 462 463 464 465 /*MC 466 SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 467 468 Level: beginner 469 470 Options Database: 471 + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter. 472 . -snes_linesearch_type <cp,l2,basic> - Line search type. 473 - -snes_ncg_monitor - Print relevant information about the ncg iteration. 474 475 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 476 gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 477 chooses the initial search direction as F(x) for the initial guess x. 478 479 480 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN 481 M*/ 482 EXTERN_C_BEGIN 483 #undef __FUNCT__ 484 #define __FUNCT__ "SNESCreate_NCG" 485 PetscErrorCode SNESCreate_NCG(SNES snes) 486 { 487 PetscErrorCode ierr; 488 SNES_NCG * neP; 489 490 PetscFunctionBegin; 491 snes->ops->destroy = SNESDestroy_NCG; 492 snes->ops->setup = SNESSetUp_NCG; 493 snes->ops->setfromoptions = SNESSetFromOptions_NCG; 494 snes->ops->view = SNESView_NCG; 495 snes->ops->solve = SNESSolve_NCG; 496 snes->ops->reset = SNESReset_NCG; 497 498 snes->usesksp = PETSC_FALSE; 499 snes->usespc = PETSC_TRUE; 500 501 if (!snes->tolerancesset) { 502 snes->max_funcs = 30000; 503 snes->max_its = 10000; 504 snes->stol = 1e-20; 505 } 506 507 ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr); 508 snes->data = (void*) neP; 509 neP->monitor = PETSC_NULL; 510 neP->type = SNES_NCG_PRP; 511 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNCGSetType_C","SNESNCGSetType_NCG", SNESNCGSetType_NCG);CHKERRQ(ierr); 512 PetscFunctionReturn(0); 513 } 514 EXTERN_C_END 515