10a844d1aSPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/ 26a6fc655SJed Brown const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0}; 3fef7b6d8SPeter Brune 4fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes) 5fef7b6d8SPeter Brune { 6fef7b6d8SPeter Brune PetscFunctionBegin; 7fef7b6d8SPeter Brune PetscFunctionReturn(0); 8fef7b6d8SPeter Brune } 9fef7b6d8SPeter Brune 109105210eSPeter Brune #define SNESLINESEARCHNCGLINEAR "ncglinear" 11bbd5d0b3SPeter Brune 12fef7b6d8SPeter Brune /* 13fef7b6d8SPeter Brune SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 14fef7b6d8SPeter Brune 15fef7b6d8SPeter Brune Input Parameter: 16fef7b6d8SPeter Brune . snes - the SNES context 17fef7b6d8SPeter Brune 18fef7b6d8SPeter Brune Application Interface Routine: SNESDestroy() 19fef7b6d8SPeter Brune */ 20fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes) 21fef7b6d8SPeter Brune { 22fef7b6d8SPeter Brune PetscErrorCode ierr; 23fef7b6d8SPeter Brune 24fef7b6d8SPeter Brune PetscFunctionBegin; 25fef7b6d8SPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 26fef7b6d8SPeter Brune PetscFunctionReturn(0); 27fef7b6d8SPeter Brune } 28fef7b6d8SPeter Brune 29fef7b6d8SPeter Brune /* 30fef7b6d8SPeter Brune SNESSetUp_NCG - Sets up the internal data structures for the later use 31fef7b6d8SPeter Brune of the SNESNCG nonlinear solver. 32fef7b6d8SPeter Brune 33fef7b6d8SPeter Brune Input Parameters: 34fef7b6d8SPeter Brune + snes - the SNES context 35fef7b6d8SPeter Brune - x - the solution vector 36fef7b6d8SPeter Brune 37fef7b6d8SPeter Brune Application Interface Routine: SNESSetUp() 38fef7b6d8SPeter Brune */ 39bbd5d0b3SPeter Brune 408cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch); 41bbd5d0b3SPeter Brune 42fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes) 43fef7b6d8SPeter Brune { 44fef7b6d8SPeter Brune PetscErrorCode ierr; 45fef7b6d8SPeter Brune 46fef7b6d8SPeter Brune PetscFunctionBegin; 47fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr); 48*efd4aadfSBarry Smith if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 496c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 50fef7b6d8SPeter Brune PetscFunctionReturn(0); 51fef7b6d8SPeter Brune } 52fef7b6d8SPeter Brune /* 5305b53524SPeter Brune SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 54fef7b6d8SPeter Brune 55fef7b6d8SPeter Brune Input Parameter: 56fef7b6d8SPeter Brune . snes - the SNES context 57fef7b6d8SPeter Brune 58fef7b6d8SPeter Brune Application Interface Routine: SNESSetFromOptions() 59fef7b6d8SPeter Brune */ 604416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes) 61fef7b6d8SPeter Brune { 62dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 63fef7b6d8SPeter Brune PetscErrorCode ierr; 6494ae4db5SBarry Smith PetscBool debug = PETSC_FALSE; 65f1c6b773SPeter Brune SNESLineSearch linesearch; 66a11d90f7SPeter Brune SNESNCGType ncgtype=ncg->type; 67fef7b6d8SPeter Brune 68fef7b6d8SPeter Brune PetscFunctionBegin; 69e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES NCG options");CHKERRQ(ierr); 700298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr); 71dfb256c7SPeter Brune if (debug) { 72ce94432eSBarry Smith ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 73dfb256c7SPeter Brune } 7494ae4db5SBarry Smith ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr); 7594ae4db5SBarry Smith ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 76fef7b6d8SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 779e764e56SPeter Brune if (!snes->linesearch) { 787601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 79*efd4aadfSBarry Smith if (!snes->npc) { 801a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 819e764e56SPeter Brune } else { 821a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 839e764e56SPeter Brune } 849e764e56SPeter Brune } 859105210eSPeter Brune ierr = SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr); 86fef7b6d8SPeter Brune PetscFunctionReturn(0); 87fef7b6d8SPeter Brune } 88fef7b6d8SPeter Brune 89fef7b6d8SPeter Brune /* 90fef7b6d8SPeter Brune SNESView_NCG - Prints info from the SNESNCG data structure. 91fef7b6d8SPeter Brune 92fef7b6d8SPeter Brune Input Parameters: 93fef7b6d8SPeter Brune + SNES - the SNES context 94fef7b6d8SPeter Brune - viewer - visualization context 95fef7b6d8SPeter Brune 96fef7b6d8SPeter Brune Application Interface Routine: SNESView() 97fef7b6d8SPeter Brune */ 98fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 99fef7b6d8SPeter Brune { 100cc3c3c0fSMatthew G. Knepley SNES_NCG *ncg = (SNES_NCG *) snes->data; 101fef7b6d8SPeter Brune PetscBool iascii; 102fef7b6d8SPeter Brune PetscErrorCode ierr; 103fef7b6d8SPeter Brune 104fef7b6d8SPeter Brune PetscFunctionBegin; 105251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 106fef7b6d8SPeter Brune if (iascii) { 107*efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type]);CHKERRQ(ierr); 108fef7b6d8SPeter Brune } 109fef7b6d8SPeter Brune PetscFunctionReturn(0); 110fef7b6d8SPeter Brune } 111fef7b6d8SPeter Brune 112f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 113ea630c6eSPeter Brune { 114ea630c6eSPeter Brune PetscScalar alpha, ptAp; 115bbd5d0b3SPeter Brune Vec X, Y, F, W; 116bbd5d0b3SPeter Brune SNES snes; 117ea630c6eSPeter Brune PetscErrorCode ierr; 118bbd5d0b3SPeter Brune PetscReal *fnorm, *xnorm, *ynorm; 119ea630c6eSPeter Brune 120ea630c6eSPeter Brune PetscFunctionBegin; 121f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 122bbd5d0b3SPeter Brune X = linesearch->vec_sol; 123bbd5d0b3SPeter Brune W = linesearch->vec_sol_new; 124bbd5d0b3SPeter Brune F = linesearch->vec_func; 125bbd5d0b3SPeter Brune Y = linesearch->vec_update; 126bbd5d0b3SPeter Brune fnorm = &linesearch->fnorm; 127bbd5d0b3SPeter Brune xnorm = &linesearch->xnorm; 128bbd5d0b3SPeter Brune ynorm = &linesearch->ynorm; 129bbd5d0b3SPeter Brune 1309105210eSPeter Brune if (!snes->jacobian) { 1319105210eSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 1329105210eSPeter Brune } 1339105210eSPeter Brune 134ea630c6eSPeter Brune /* 135ea630c6eSPeter Brune 136d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 137ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 138ea630c6eSPeter Brune */ 139d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);CHKERRQ(ierr); 140ea630c6eSPeter Brune ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 141ea630c6eSPeter Brune ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 142ea630c6eSPeter Brune ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 143ea630c6eSPeter Brune alpha = alpha / ptAp; 1449105210eSPeter Brune ierr = VecAXPY(X,-alpha,Y);CHKERRQ(ierr); 145bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 146bbd5d0b3SPeter Brune 147bbd5d0b3SPeter Brune ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 148bbd5d0b3SPeter Brune ierr = VecNorm(X,NORM_2,xnorm);CHKERRQ(ierr); 149bbd5d0b3SPeter Brune ierr = VecNorm(Y,NORM_2,ynorm);CHKERRQ(ierr); 150ea630c6eSPeter Brune PetscFunctionReturn(0); 151ea630c6eSPeter Brune } 152ea630c6eSPeter Brune 1538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 154bbd5d0b3SPeter Brune { 155bbd5d0b3SPeter Brune PetscFunctionBegin; 156f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 1570298fd71SBarry Smith linesearch->ops->destroy = NULL; 1580298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1590298fd71SBarry Smith linesearch->ops->reset = NULL; 1600298fd71SBarry Smith linesearch->ops->view = NULL; 1610298fd71SBarry Smith linesearch->ops->setup = NULL; 162bbd5d0b3SPeter Brune PetscFunctionReturn(0); 163bbd5d0b3SPeter Brune } 164bbd5d0b3SPeter Brune 16588d374b2SPeter Brune /* 16688d374b2SPeter Brune 16788d374b2SPeter Brune Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 16888d374b2SPeter Brune 16988d374b2SPeter Brune */ 17067392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) 17167392de3SBarry Smith { 17288d374b2SPeter Brune PetscErrorCode ierr; 17388d374b2SPeter Brune PetscScalar ftf, ftg, fty, h; 17467392de3SBarry Smith 17588d374b2SPeter Brune PetscFunctionBegin; 17688d374b2SPeter Brune ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 17788d374b2SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 17888d374b2SPeter Brune h = 1e-5*fty / fty; 17988d374b2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 18088d374b2SPeter Brune ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 18188d374b2SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 18288d374b2SPeter Brune ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 18388d374b2SPeter Brune *ytJtf = (ftg - ftf) / h; 18488d374b2SPeter Brune PetscFunctionReturn(0); 18588d374b2SPeter Brune } 18688d374b2SPeter Brune 1870a844d1aSPeter Brune /*@ 1880a844d1aSPeter Brune SNESNCGSetType - Sets the conjugate update type for SNESNCG. 1890a844d1aSPeter Brune 1900a844d1aSPeter Brune Logically Collective on SNES 1910a844d1aSPeter Brune 1920a844d1aSPeter Brune Input Parameters: 1930a844d1aSPeter Brune + snes - the iterative context 1940a844d1aSPeter Brune - btype - update type 1950a844d1aSPeter Brune 1960a844d1aSPeter Brune Options Database: 1970a844d1aSPeter Brune . -snes_ncg_type<prp,fr,hs,dy,cd> 1980a844d1aSPeter Brune 1990a844d1aSPeter Brune Level: intermediate 2000a844d1aSPeter Brune 2010a844d1aSPeter Brune SNESNCGSelectTypes: 2020a844d1aSPeter Brune + SNES_NCG_FR - Fletcher-Reeves update 2030a844d1aSPeter Brune . SNES_NCG_PRP - Polak-Ribiere-Polyak update 2040a844d1aSPeter Brune . SNES_NCG_HS - Hestenes-Steifel update 2050a844d1aSPeter Brune . SNES_NCG_DY - Dai-Yuan update 2060a844d1aSPeter Brune - SNES_NCG_CD - Conjugate Descent update 2070a844d1aSPeter Brune 2080a844d1aSPeter Brune Notes: 2090a844d1aSPeter Brune PRP is the default, and the only one that tolerates generalized search directions. 2100a844d1aSPeter Brune 2110a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set 2120a844d1aSPeter Brune @*/ 2130adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 2140adebc6cSBarry Smith { 2150a844d1aSPeter Brune PetscErrorCode ierr; 2166e111a19SKarl Rupp 2170a844d1aSPeter Brune PetscFunctionBegin; 2180a844d1aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2190a844d1aSPeter Brune ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr); 2200a844d1aSPeter Brune PetscFunctionReturn(0); 2210a844d1aSPeter Brune } 2220a844d1aSPeter Brune 2230adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 2240adebc6cSBarry Smith { 2250a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 2266e111a19SKarl Rupp 2270a844d1aSPeter Brune PetscFunctionBegin; 2280a844d1aSPeter Brune ncg->type = btype; 2290a844d1aSPeter Brune PetscFunctionReturn(0); 2300a844d1aSPeter Brune } 2310a844d1aSPeter Brune 232fef7b6d8SPeter Brune /* 233fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 234fef7b6d8SPeter Brune 235fef7b6d8SPeter Brune Input Parameters: 236fef7b6d8SPeter Brune . snes - the SNES context 237fef7b6d8SPeter Brune 238fef7b6d8SPeter Brune Output Parameter: 239fef7b6d8SPeter Brune . outits - number of iterations until termination 240fef7b6d8SPeter Brune 241fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 242fef7b6d8SPeter Brune */ 243fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes) 244fef7b6d8SPeter Brune { 245dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 24632cc618eSPeter Brune Vec X,dX,lX,F,dXold; 247bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0; 24832cc618eSPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 249fef7b6d8SPeter Brune PetscInt maxits, i; 250fef7b6d8SPeter Brune PetscErrorCode ierr; 251422a814eSBarry Smith SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 252f1c6b773SPeter Brune SNESLineSearch linesearch; 253b7281c8aSPeter Brune SNESConvergedReason reason; 25488d374b2SPeter Brune 255fef7b6d8SPeter Brune PetscFunctionBegin; 2566c4ed002SBarry Smith 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); 257c579b300SPatrick Farrell 258fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 259fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 260fef7b6d8SPeter Brune 261fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 262fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 26332cc618eSPeter Brune dXold = snes->work[0]; /* The previous iterate of X */ 26488d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */ 265169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 266fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 267fef7b6d8SPeter Brune 2687601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 2699e764e56SPeter Brune 270e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 271fef7b6d8SPeter Brune snes->iter = 0; 272fef7b6d8SPeter Brune snes->norm = 0.; 273e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 274fef7b6d8SPeter Brune 275bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */ 276a71552e2SPeter Brune 277*efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 278be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 279*efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 280b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 281b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 282b7281c8aSPeter Brune PetscFunctionReturn(0); 283b7281c8aSPeter Brune } 284a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 285a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 286a71552e2SPeter Brune } else { 287e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 288fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2891aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 290c1c75074SPeter Brune 291fef7b6d8SPeter Brune /* convergence test */ 292a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 293422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 294a71552e2SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 295a71552e2SPeter Brune } 296*efd4aadfSBarry Smith if (snes->npc) { 297a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 298be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 299*efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 300b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 301b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 302b7281c8aSPeter Brune PetscFunctionReturn(0); 303b7281c8aSPeter Brune } 304a71552e2SPeter Brune } 305a71552e2SPeter Brune } 306a71552e2SPeter Brune ierr = VecCopy(dX,lX);CHKERRQ(ierr); 30732cc618eSPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 308a71552e2SPeter Brune 309a71552e2SPeter Brune 310e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 311fef7b6d8SPeter Brune snes->norm = fnorm; 312e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 313a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 314fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 315fef7b6d8SPeter Brune 316fef7b6d8SPeter Brune /* test convergence */ 317fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 318fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 319fef7b6d8SPeter Brune 320fef7b6d8SPeter Brune /* Call general purpose update function */ 321fef7b6d8SPeter Brune if (snes->ops->update) { 322fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 323fef7b6d8SPeter Brune } 324fef7b6d8SPeter Brune 325fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 326bbd5d0b3SPeter Brune 32709c08436SPeter Brune for (i = 1; i < maxits + 1; i++) { 32829c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 3290a844d1aSPeter Brune if (ncg->type != SNES_NCG_FR) { 33032cc618eSPeter Brune ierr = VecCopy(dX, dXold);CHKERRQ(ierr); 33188d374b2SPeter Brune } 332f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr); 333422a814eSBarry Smith ierr = SNESLineSearchGetReason(linesearch, &lsresult);CHKERRQ(ierr); 334422a814eSBarry Smith ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 335422a814eSBarry Smith if (lsresult) { 336fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 337fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 3385d10c001SPeter Brune PetscFunctionReturn(0); 339fef7b6d8SPeter Brune } 340fef7b6d8SPeter Brune } 341fef7b6d8SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 342fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 3435d10c001SPeter Brune PetscFunctionReturn(0); 344fef7b6d8SPeter Brune } 345fef7b6d8SPeter Brune /* Monitor convergence */ 346e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 347169e2be8SPeter Brune snes->iter = i; 348fef7b6d8SPeter Brune snes->norm = fnorm; 349e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 350a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 351fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 352302dbbaeSPeter Brune 353fef7b6d8SPeter Brune /* Test for convergence */ 354d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3555d10c001SPeter Brune if (snes->reason) PetscFunctionReturn(0); 356302dbbaeSPeter Brune 357302dbbaeSPeter Brune /* Call general purpose update function */ 358302dbbaeSPeter Brune if (snes->ops->update) { 359302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 360302dbbaeSPeter Brune } 361*efd4aadfSBarry Smith if (snes->npc) { 362a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 363be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 364*efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 365b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 366b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 367b7281c8aSPeter Brune PetscFunctionReturn(0); 368b7281c8aSPeter Brune } 369a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 370a71552e2SPeter Brune } else { 371be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 372*efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 373b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 374b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 375b7281c8aSPeter Brune PetscFunctionReturn(0); 376b7281c8aSPeter Brune } 377302dbbaeSPeter Brune } 378a3685310SPeter Brune } else { 379a3685310SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 380302dbbaeSPeter Brune } 381302dbbaeSPeter Brune 38229c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 3830a844d1aSPeter Brune switch (ncg->type) { 3840a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */ 38532cc618eSPeter Brune dXolddotdXold= dXdotdX; 38632cc618eSPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 38732cc618eSPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold); 388dfb256c7SPeter Brune break; 3890a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 39032cc618eSPeter Brune dXolddotdXold= dXdotdX; 39132cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr); 39232cc618eSPeter Brune ierr = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 39332cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 39432cc618eSPeter Brune ierr = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 39532cc618eSPeter Brune beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 396dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 397dfb256c7SPeter Brune break; 3980a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */ 39932cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 40032cc618eSPeter Brune ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 40132cc618eSPeter Brune ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 40232cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 40332cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 40432cc618eSPeter Brune ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 40532cc618eSPeter Brune ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 40632cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 40732cc618eSPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 408dfb256c7SPeter Brune break; 4090a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */ 41032cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 41132cc618eSPeter Brune ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 41232cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 41332cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 41432cc618eSPeter Brune ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 41532cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 41632cc618eSPeter Brune beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr); 417dfb256c7SPeter Brune break; 4180a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */ 41932cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 42032cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 42132cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 42232cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 42332cc618eSPeter Brune beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr); 424dfb256c7SPeter Brune break; 425dfb256c7SPeter Brune } 426dfb256c7SPeter Brune if (ncg->monitor) { 4277904a332SBarry Smith ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);CHKERRQ(ierr); 428dfb256c7SPeter Brune } 429dfb256c7SPeter Brune ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 430fef7b6d8SPeter Brune } 431fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 432fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 433fef7b6d8SPeter Brune PetscFunctionReturn(0); 434fef7b6d8SPeter Brune } 435fef7b6d8SPeter Brune 4360a844d1aSPeter Brune 4370a844d1aSPeter Brune 438fef7b6d8SPeter Brune /*MC 439fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 440fef7b6d8SPeter Brune 441fef7b6d8SPeter Brune Level: beginner 442fef7b6d8SPeter Brune 443fef7b6d8SPeter Brune Options Database: 4444f02bc6aSBarry Smith + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp. 44541a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type. 4461867fe5bSPeter Brune - -snes_ncg_monitor - Print relevant information about the ncg iteration. 447fef7b6d8SPeter Brune 448fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 449fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 450fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 451fef7b6d8SPeter Brune 4522d547940SBarry Smith Only supports left non-linear preconditioning. 4532d547940SBarry Smith 4544f02bc6aSBarry Smith References: 45596a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 4564f02bc6aSBarry Smith SIAM Review, 57(4), 2015 457cc3c3c0fSMatthew G. Knepley 458fef7b6d8SPeter Brune 45904d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN 460fef7b6d8SPeter Brune M*/ 4618cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 462fef7b6d8SPeter Brune { 463fef7b6d8SPeter Brune PetscErrorCode ierr; 464ea630c6eSPeter Brune SNES_NCG * neP; 465fef7b6d8SPeter Brune 466fef7b6d8SPeter Brune PetscFunctionBegin; 467fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 468fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 469fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 470fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 471fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 472fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 473fef7b6d8SPeter Brune 474fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 475*efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 476*efd4aadfSBarry Smith snes->npcside = PC_LEFT; 477fef7b6d8SPeter Brune 4784fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 4794fc747eaSLawrence Mitchell 48088976e71SPeter Brune if (!snes->tolerancesset) { 4810e444f03SPeter Brune snes->max_funcs = 30000; 4820e444f03SPeter Brune snes->max_its = 10000; 483c522fa08SPeter Brune snes->stol = 1e-20; 48488976e71SPeter Brune } 4850e444f03SPeter Brune 486b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 487fef7b6d8SPeter Brune snes->data = (void*) neP; 4880298fd71SBarry Smith neP->monitor = NULL; 4890a844d1aSPeter Brune neP->type = SNES_NCG_PRP; 490bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr); 491fef7b6d8SPeter Brune PetscFunctionReturn(0); 492fef7b6d8SPeter Brune } 493