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 static PetscErrorCode SNESReset_NCG(SNES snes) 5 { 6 PetscFunctionBegin; 7 PetscFunctionReturn(0); 8 } 9 10 /* 11 SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 12 13 Input Parameter: 14 . snes - the SNES context 15 16 Application Interface Routine: SNESDestroy() 17 */ 18 static PetscErrorCode SNESDestroy_NCG(SNES snes) 19 { 20 PetscErrorCode ierr; 21 22 PetscFunctionBegin; 23 ierr = PetscFree(snes->data);CHKERRQ(ierr); 24 PetscFunctionReturn(0); 25 } 26 27 /* 28 SNESSetUp_NCG - Sets up the internal data structures for the later use 29 of the SNESNCG nonlinear solver. 30 31 Input Parameters: 32 + snes - the SNES context 33 - x - the solution vector 34 35 Application Interface Routine: SNESSetUp() 36 */ 37 38 static PetscErrorCode SNESSetUp_NCG(SNES snes) 39 { 40 PetscErrorCode ierr; 41 42 PetscFunctionBegin; 43 ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr); 44 if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 45 if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 46 PetscFunctionReturn(0); 47 } 48 49 static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 50 { 51 PetscScalar alpha, ptAp; 52 Vec X, Y, F, W; 53 SNES snes; 54 PetscErrorCode ierr; 55 PetscReal *fnorm, *xnorm, *ynorm; 56 57 PetscFunctionBegin; 58 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 59 X = linesearch->vec_sol; 60 W = linesearch->vec_sol_new; 61 F = linesearch->vec_func; 62 Y = linesearch->vec_update; 63 fnorm = &linesearch->fnorm; 64 xnorm = &linesearch->xnorm; 65 ynorm = &linesearch->ynorm; 66 67 if (!snes->jacobian) { 68 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 69 } 70 71 /* 72 73 The exact step size for unpreconditioned linear CG is just: 74 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 75 */ 76 ierr = SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);CHKERRQ(ierr); 77 ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 78 ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 79 ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 80 alpha = alpha / ptAp; 81 ierr = VecAXPY(X,-alpha,Y);CHKERRQ(ierr); 82 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 83 84 ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 85 ierr = VecNorm(X,NORM_2,xnorm);CHKERRQ(ierr); 86 ierr = VecNorm(Y,NORM_2,ynorm);CHKERRQ(ierr); 87 PetscFunctionReturn(0); 88 } 89 90 /*MC 91 SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG 92 93 This line search uses the length "as if" the problem is linear (that is what is computed by the linear CG method) using the Jacobian of the function. 94 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) where r (f) is the current residual (function value), p (y) is the current search direction. 95 96 Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian 97 98 This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone. 99 100 Level: advanced 101 102 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 103 M*/ 104 105 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 106 { 107 PetscFunctionBegin; 108 linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 109 linesearch->ops->destroy = NULL; 110 linesearch->ops->setfromoptions = NULL; 111 linesearch->ops->reset = NULL; 112 linesearch->ops->view = NULL; 113 linesearch->ops->setup = NULL; 114 PetscFunctionReturn(0); 115 } 116 117 /* 118 SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 119 120 Input Parameter: 121 . snes - the SNES context 122 123 Application Interface Routine: SNESSetFromOptions() 124 */ 125 static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes) 126 { 127 SNES_NCG *ncg = (SNES_NCG*)snes->data; 128 PetscErrorCode ierr; 129 PetscBool debug = PETSC_FALSE; 130 SNESNCGType ncgtype=ncg->type; 131 SNESLineSearch linesearch; 132 133 PetscFunctionBegin; 134 ierr = PetscOptionsHead(PetscOptionsObject,"SNES NCG options");CHKERRQ(ierr); 135 ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr); 136 if (debug) { 137 ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 138 } 139 ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr); 140 ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 141 ierr = PetscOptionsTail();CHKERRQ(ierr); 142 if (!snes->linesearch) { 143 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 144 if (!snes->npc) { 145 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 146 } else { 147 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 148 } 149 } 150 PetscFunctionReturn(0); 151 } 152 153 /* 154 SNESView_NCG - Prints info from the SNESNCG data structure. 155 156 Input Parameters: 157 + SNES - the SNES context 158 - viewer - visualization context 159 160 Application Interface Routine: SNESView() 161 */ 162 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 163 { 164 SNES_NCG *ncg = (SNES_NCG *) snes->data; 165 PetscBool iascii; 166 PetscErrorCode ierr; 167 168 PetscFunctionBegin; 169 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 170 if (iascii) { 171 ierr = PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type]);CHKERRQ(ierr); 172 } 173 PetscFunctionReturn(0); 174 } 175 176 /* 177 178 Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 179 180 This routine is not currently used. I don't know what its intended purpose is. 181 182 Note it has a hardwired differencing parameter of 1e-5 183 184 */ 185 PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) 186 { 187 PetscErrorCode ierr; 188 PetscScalar ftf, ftg, fty, h; 189 190 PetscFunctionBegin; 191 ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 192 ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 193 h = 1e-5*fty / fty; 194 ierr = VecCopy(X, W);CHKERRQ(ierr); 195 ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 196 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 197 ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 198 *ytJtf = (ftg - ftf) / h; 199 PetscFunctionReturn(0); 200 } 201 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> -strategy for selecting algorithm for computing beta 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 SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions. 225 226 It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner, 227 that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()? 228 229 @*/ 230 PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 231 { 232 PetscErrorCode ierr; 233 234 PetscFunctionBegin; 235 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 236 ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr); 237 PetscFunctionReturn(0); 238 } 239 240 static 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 static PetscErrorCode SNESSolve_NCG(SNES snes) 261 { 262 SNES_NCG *ncg = (SNES_NCG*)snes->data; 263 Vec X,dX,lX,F,dXold; 264 PetscReal fnorm, ynorm, xnorm, beta = 0.0; 265 PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 266 PetscInt maxits, i; 267 PetscErrorCode ierr; 268 SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 269 SNESLineSearch linesearch; 270 SNESConvergedReason reason; 271 272 PetscFunctionBegin; 273 if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 274 275 ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 276 snes->reason = SNES_CONVERGED_ITERATING; 277 278 maxits = snes->max_its; /* maximum number of iterations */ 279 X = snes->vec_sol; /* X^n */ 280 dXold = snes->work[0]; /* The previous iterate of X */ 281 dX = snes->work[1]; /* the preconditioned direction */ 282 lX = snes->vec_sol_update; /* the conjugate direction */ 283 F = snes->vec_func; /* residual vector */ 284 285 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 286 287 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 288 snes->iter = 0; 289 snes->norm = 0.; 290 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 291 292 /* compute the initial function and preconditioned update dX */ 293 294 if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 295 ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 296 ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 297 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 298 snes->reason = SNES_DIVERGED_INNER; 299 PetscFunctionReturn(0); 300 } 301 ierr = VecCopy(dX,F);CHKERRQ(ierr); 302 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 303 } else { 304 if (!snes->vec_func_init_set) { 305 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 306 } else snes->vec_func_init_set = PETSC_FALSE; 307 308 /* convergence test */ 309 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 310 SNESCheckFunctionNorm(snes,fnorm); 311 ierr = VecCopy(F,dX);CHKERRQ(ierr); 312 } 313 if (snes->npc) { 314 if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 315 ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 316 ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 317 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 318 snes->reason = SNES_DIVERGED_INNER; 319 PetscFunctionReturn(0); 320 } 321 } 322 } 323 ierr = VecCopy(dX,lX);CHKERRQ(ierr); 324 ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 325 326 327 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 328 snes->norm = fnorm; 329 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 330 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 331 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 332 333 /* test convergence */ 334 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 335 if (snes->reason) PetscFunctionReturn(0); 336 337 /* Call general purpose update function */ 338 if (snes->ops->update) { 339 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 340 } 341 342 /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 343 344 for (i = 1; i < maxits + 1; i++) { 345 /* some update types require the old update direction or conjugate direction */ 346 if (ncg->type != SNES_NCG_FR) { 347 ierr = VecCopy(dX, dXold);CHKERRQ(ierr); 348 } 349 ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr); 350 ierr = SNESLineSearchGetReason(linesearch, &lsresult);CHKERRQ(ierr); 351 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 352 if (lsresult) { 353 if (++snes->numFailures >= snes->maxFailures) { 354 snes->reason = SNES_DIVERGED_LINE_SEARCH; 355 PetscFunctionReturn(0); 356 } 357 } 358 if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 359 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 360 PetscFunctionReturn(0); 361 } 362 /* Monitor convergence */ 363 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 364 snes->iter = i; 365 snes->norm = fnorm; 366 snes->xnorm = xnorm; 367 snes->ynorm = ynorm; 368 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 369 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 370 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 371 372 /* Test for convergence */ 373 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 374 if (snes->reason) PetscFunctionReturn(0); 375 376 /* Call general purpose update function */ 377 if (snes->ops->update) { 378 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 379 } 380 if (snes->npc) { 381 if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 382 ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 383 ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 384 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 385 snes->reason = SNES_DIVERGED_INNER; 386 PetscFunctionReturn(0); 387 } 388 ierr = VecCopy(dX,F);CHKERRQ(ierr); 389 } else { 390 ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 391 ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 392 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 393 snes->reason = SNES_DIVERGED_INNER; 394 PetscFunctionReturn(0); 395 } 396 } 397 } else { 398 ierr = VecCopy(F,dX);CHKERRQ(ierr); 399 } 400 401 /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 402 switch (ncg->type) { 403 case SNES_NCG_FR: /* Fletcher-Reeves */ 404 dXolddotdXold= dXdotdX; 405 ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 406 beta = PetscRealPart(dXdotdX / dXolddotdXold); 407 break; 408 case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 409 dXolddotdXold= dXdotdX; 410 ierr = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr); 411 ierr = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 412 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 413 ierr = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 414 beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 415 if (beta < 0.0) beta = 0.0; /* restart */ 416 break; 417 case SNES_NCG_HS: /* Hestenes-Stiefel */ 418 ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 419 ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 420 ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 421 ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 422 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 423 ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 424 ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 425 ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 426 beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 427 break; 428 case SNES_NCG_DY: /* Dai-Yuan */ 429 ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 430 ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 431 ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 432 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 433 ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 434 ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 435 beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr); 436 break; 437 case SNES_NCG_CD: /* Conjugate Descent */ 438 ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 439 ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 440 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 441 ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 442 beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr); 443 break; 444 } 445 if (ncg->monitor) { 446 ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);CHKERRQ(ierr); 447 } 448 ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 449 } 450 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 451 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 452 PetscFunctionReturn(0); 453 } 454 455 /*MC 456 SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 457 458 Level: beginner 459 460 Options Database: 461 + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp. 462 . -snes_linesearch_type <cp,l2,basic> - Line search type. 463 - -snes_ncg_monitor - Print the beta values used in the ncg iteration, . 464 465 Notes: 466 This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 467 gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 468 chooses the initial search direction as F(x) for the initial guess x. 469 470 Only supports left non-linear preconditioning. 471 472 Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then 473 SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR 474 475 References: 476 . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 477 SIAM Review, 57(4), 2015 478 479 480 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType() 481 M*/ 482 PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 483 { 484 PetscErrorCode ierr; 485 SNES_NCG * neP; 486 487 PetscFunctionBegin; 488 snes->ops->destroy = SNESDestroy_NCG; 489 snes->ops->setup = SNESSetUp_NCG; 490 snes->ops->setfromoptions = SNESSetFromOptions_NCG; 491 snes->ops->view = SNESView_NCG; 492 snes->ops->solve = SNESSolve_NCG; 493 snes->ops->reset = SNESReset_NCG; 494 495 snes->usesksp = PETSC_FALSE; 496 snes->usesnpc = PETSC_TRUE; 497 snes->npcside = PC_LEFT; 498 499 snes->alwayscomputesfinalresidual = 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,&neP);CHKERRQ(ierr); 508 snes->data = (void*) neP; 509 neP->monitor = NULL; 510 neP->type = SNES_NCG_PRP; 511 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr); 512 PetscFunctionReturn(0); 513 } 514