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