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