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 if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 55 if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 56 PetscFunctionReturn(0); 57 } 58 /* 59 SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 60 61 Input Parameter: 62 . snes - the SNES context 63 64 Application Interface Routine: SNESSetFromOptions() 65 */ 66 #undef __FUNCT__ 67 #define __FUNCT__ "SNESSetFromOptions_NCG" 68 static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes) 69 { 70 SNES_NCG *ncg = (SNES_NCG*)snes->data; 71 PetscErrorCode ierr; 72 PetscBool debug = PETSC_FALSE; 73 SNESLineSearch linesearch; 74 SNESNCGType ncgtype=ncg->type; 75 76 PetscFunctionBegin; 77 ierr = PetscOptionsHead(PetscOptionsObject,"SNES NCG options");CHKERRQ(ierr); 78 ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr); 79 if (debug) { 80 ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 81 } 82 ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr); 83 ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 84 ierr = PetscOptionsTail();CHKERRQ(ierr); 85 if (!snes->linesearch) { 86 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 87 if (!snes->pc) { 88 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 89 } else { 90 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 91 } 92 } 93 ierr = SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr); 94 PetscFunctionReturn(0); 95 } 96 97 /* 98 SNESView_NCG - Prints info from the SNESNCG data structure. 99 100 Input Parameters: 101 + SNES - the SNES context 102 - viewer - visualization context 103 104 Application Interface Routine: SNESView() 105 */ 106 #undef __FUNCT__ 107 #define __FUNCT__ "SNESView_NCG" 108 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 109 { 110 SNES_NCG *ncg = (SNES_NCG *) snes->data; 111 PetscBool iascii; 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 116 if (iascii) { 117 ierr = PetscViewerASCIIPrintf(viewer, "ncg type: %s\n", SNESNCGTypes[ncg->type]);CHKERRQ(ierr); 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 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 if (!snes->jacobian) { 143 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 144 } 145 146 /* 147 148 The exact step size for unpreconditioned linear CG is just: 149 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 150 */ 151 ierr = SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);CHKERRQ(ierr); 152 ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 153 ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 154 ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 155 alpha = alpha / ptAp; 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 #undef __FUNCT__ 166 #define __FUNCT__ "SNESLineSearchCreate_NCGLinear" 167 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 168 { 169 PetscFunctionBegin; 170 linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 171 linesearch->ops->destroy = NULL; 172 linesearch->ops->setfromoptions = NULL; 173 linesearch->ops->reset = NULL; 174 linesearch->ops->view = NULL; 175 linesearch->ops->setup = NULL; 176 PetscFunctionReturn(0); 177 } 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 #undef __FUNCT__ 242 #define __FUNCT__ "SNESNCGSetType_NCG" 243 PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 244 { 245 SNES_NCG *ncg = (SNES_NCG*)snes->data; 246 247 PetscFunctionBegin; 248 ncg->type = btype; 249 PetscFunctionReturn(0); 250 } 251 252 /* 253 SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 254 255 Input Parameters: 256 . snes - the SNES context 257 258 Output Parameter: 259 . outits - number of iterations until termination 260 261 Application Interface Routine: SNESSolve() 262 */ 263 #undef __FUNCT__ 264 #define __FUNCT__ "SNESSolve_NCG" 265 PetscErrorCode SNESSolve_NCG(SNES snes) 266 { 267 SNES_NCG *ncg = (SNES_NCG*)snes->data; 268 Vec X,dX,lX,F,dXold; 269 PetscReal fnorm, ynorm, xnorm, beta = 0.0; 270 PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 271 PetscInt maxits, i; 272 PetscErrorCode ierr; 273 SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 274 SNESLineSearch linesearch; 275 SNESConvergedReason reason; 276 277 PetscFunctionBegin; 278 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); 279 280 ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 281 snes->reason = SNES_CONVERGED_ITERATING; 282 283 maxits = snes->max_its; /* maximum number of iterations */ 284 X = snes->vec_sol; /* X^n */ 285 dXold = 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 290 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 291 292 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 293 snes->iter = 0; 294 snes->norm = 0.; 295 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 296 297 /* compute the initial function and preconditioned update dX */ 298 299 if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 300 ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 301 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 302 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 303 snes->reason = SNES_DIVERGED_INNER; 304 PetscFunctionReturn(0); 305 } 306 ierr = VecCopy(dX,F);CHKERRQ(ierr); 307 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 308 } else { 309 if (!snes->vec_func_init_set) { 310 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 311 } else snes->vec_func_init_set = PETSC_FALSE; 312 313 /* convergence test */ 314 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 315 SNESCheckFunctionNorm(snes,fnorm); 316 ierr = VecCopy(F,dX);CHKERRQ(ierr); 317 } 318 if (snes->pc) { 319 if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 320 ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 321 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 322 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 323 snes->reason = SNES_DIVERGED_INNER; 324 PetscFunctionReturn(0); 325 } 326 } 327 } 328 ierr = VecCopy(dX,lX);CHKERRQ(ierr); 329 ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 330 331 332 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 333 snes->norm = fnorm; 334 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 335 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 336 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 337 338 /* test convergence */ 339 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 340 if (snes->reason) PetscFunctionReturn(0); 341 342 /* Call general purpose update function */ 343 if (snes->ops->update) { 344 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 345 } 346 347 /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 348 349 for (i = 1; i < maxits + 1; i++) { 350 /* some update types require the old update direction or conjugate direction */ 351 if (ncg->type != SNES_NCG_FR) { 352 ierr = VecCopy(dX, dXold);CHKERRQ(ierr); 353 } 354 ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr); 355 ierr = SNESLineSearchGetReason(linesearch, &lsresult);CHKERRQ(ierr); 356 ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 357 if (lsresult) { 358 if (++snes->numFailures >= snes->maxFailures) { 359 snes->reason = SNES_DIVERGED_LINE_SEARCH; 360 PetscFunctionReturn(0); 361 } 362 } 363 if (snes->nfuncs >= snes->max_funcs) { 364 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 365 PetscFunctionReturn(0); 366 } 367 /* Monitor convergence */ 368 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 369 snes->iter = i; 370 snes->norm = fnorm; 371 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 372 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 373 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 374 375 /* Test for convergence */ 376 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 377 if (snes->reason) PetscFunctionReturn(0); 378 379 /* Call general purpose update function */ 380 if (snes->ops->update) { 381 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 382 } 383 if (snes->pc) { 384 if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 385 ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 386 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 387 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 388 snes->reason = SNES_DIVERGED_INNER; 389 PetscFunctionReturn(0); 390 } 391 ierr = VecCopy(dX,F);CHKERRQ(ierr); 392 } else { 393 ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 394 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 395 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 396 snes->reason = SNES_DIVERGED_INNER; 397 PetscFunctionReturn(0); 398 } 399 } 400 } else { 401 ierr = VecCopy(F,dX);CHKERRQ(ierr); 402 } 403 404 /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 405 switch (ncg->type) { 406 case SNES_NCG_FR: /* Fletcher-Reeves */ 407 dXolddotdXold= dXdotdX; 408 ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 409 beta = PetscRealPart(dXdotdX / dXolddotdXold); 410 break; 411 case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 412 dXolddotdXold= dXdotdX; 413 ierr = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr); 414 ierr = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 415 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 416 ierr = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 417 beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 418 if (beta < 0.0) beta = 0.0; /* restart */ 419 break; 420 case SNES_NCG_HS: /* Hestenes-Stiefel */ 421 ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 422 ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 423 ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 424 ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 425 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 426 ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 427 ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 428 ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 429 beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 430 break; 431 case SNES_NCG_DY: /* Dai-Yuan */ 432 ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 433 ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 434 ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 435 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 436 ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 437 ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 438 beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr); 439 break; 440 case SNES_NCG_CD: /* Conjugate Descent */ 441 ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 442 ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 443 ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 444 ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 445 beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr); 446 break; 447 } 448 if (ncg->monitor) { 449 ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);CHKERRQ(ierr); 450 } 451 ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 452 } 453 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 454 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 455 PetscFunctionReturn(0); 456 } 457 458 459 460 /*MC 461 SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 462 463 Level: beginner 464 465 Options Database: 466 + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp. 467 . -snes_linesearch_type <cp,l2,basic> - Line search type. 468 - -snes_ncg_monitor - Print relevant information about the ncg iteration. 469 470 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 471 gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 472 chooses the initial search direction as F(x) for the initial guess x. 473 474 References: 475 . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 476 SIAM Review, 57(4), 2015 477 478 479 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN 480 M*/ 481 #undef __FUNCT__ 482 #define __FUNCT__ "SNESCreate_NCG" 483 PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 484 { 485 PetscErrorCode ierr; 486 SNES_NCG * neP; 487 488 PetscFunctionBegin; 489 snes->ops->destroy = SNESDestroy_NCG; 490 snes->ops->setup = SNESSetUp_NCG; 491 snes->ops->setfromoptions = SNESSetFromOptions_NCG; 492 snes->ops->view = SNESView_NCG; 493 snes->ops->solve = SNESSolve_NCG; 494 snes->ops->reset = SNESReset_NCG; 495 496 snes->usesksp = PETSC_FALSE; 497 snes->usespc = PETSC_TRUE; 498 snes->pcside = PC_LEFT; 499 500 if (!snes->tolerancesset) { 501 snes->max_funcs = 30000; 502 snes->max_its = 10000; 503 snes->stol = 1e-20; 504 } 505 506 ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 507 snes->data = (void*) neP; 508 neP->monitor = NULL; 509 neP->type = SNES_NCG_PRP; 510 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr); 511 PetscFunctionReturn(0); 512 } 513