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