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 4*b5badacbSBarry Smith static PetscErrorCode SNESReset_NCG(SNES snes) 5fef7b6d8SPeter Brune { 6fef7b6d8SPeter Brune PetscFunctionBegin; 7fef7b6d8SPeter Brune PetscFunctionReturn(0); 8fef7b6d8SPeter Brune } 9fef7b6d8SPeter Brune 10fef7b6d8SPeter Brune /* 11fef7b6d8SPeter Brune SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 12fef7b6d8SPeter Brune 13fef7b6d8SPeter Brune Input Parameter: 14fef7b6d8SPeter Brune . snes - the SNES context 15fef7b6d8SPeter Brune 16fef7b6d8SPeter Brune Application Interface Routine: SNESDestroy() 17fef7b6d8SPeter Brune */ 18*b5badacbSBarry Smith static PetscErrorCode SNESDestroy_NCG(SNES snes) 19fef7b6d8SPeter Brune { 20fef7b6d8SPeter Brune PetscErrorCode ierr; 21fef7b6d8SPeter Brune 22fef7b6d8SPeter Brune PetscFunctionBegin; 23fef7b6d8SPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 24fef7b6d8SPeter Brune PetscFunctionReturn(0); 25fef7b6d8SPeter Brune } 26fef7b6d8SPeter Brune 27fef7b6d8SPeter Brune /* 28fef7b6d8SPeter Brune SNESSetUp_NCG - Sets up the internal data structures for the later use 29fef7b6d8SPeter Brune of the SNESNCG nonlinear solver. 30fef7b6d8SPeter Brune 31fef7b6d8SPeter Brune Input Parameters: 32fef7b6d8SPeter Brune + snes - the SNES context 33fef7b6d8SPeter Brune - x - the solution vector 34fef7b6d8SPeter Brune 35fef7b6d8SPeter Brune Application Interface Routine: SNESSetUp() 36fef7b6d8SPeter Brune */ 37bbd5d0b3SPeter Brune 38*b5badacbSBarry Smith static PetscErrorCode SNESSetUp_NCG(SNES snes) 39fef7b6d8SPeter Brune { 40fef7b6d8SPeter Brune PetscErrorCode ierr; 41fef7b6d8SPeter Brune 42fef7b6d8SPeter Brune PetscFunctionBegin; 43fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr); 44efd4aadfSBarry Smith if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 456c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 46fef7b6d8SPeter Brune PetscFunctionReturn(0); 47fef7b6d8SPeter Brune } 48fef7b6d8SPeter Brune 49*b5badacbSBarry Smith static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 50ea630c6eSPeter Brune { 51ea630c6eSPeter Brune PetscScalar alpha, ptAp; 52bbd5d0b3SPeter Brune Vec X, Y, F, W; 53bbd5d0b3SPeter Brune SNES snes; 54ea630c6eSPeter Brune PetscErrorCode ierr; 55bbd5d0b3SPeter Brune PetscReal *fnorm, *xnorm, *ynorm; 56ea630c6eSPeter Brune 57ea630c6eSPeter Brune PetscFunctionBegin; 58f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 59bbd5d0b3SPeter Brune X = linesearch->vec_sol; 60bbd5d0b3SPeter Brune W = linesearch->vec_sol_new; 61bbd5d0b3SPeter Brune F = linesearch->vec_func; 62bbd5d0b3SPeter Brune Y = linesearch->vec_update; 63bbd5d0b3SPeter Brune fnorm = &linesearch->fnorm; 64bbd5d0b3SPeter Brune xnorm = &linesearch->xnorm; 65bbd5d0b3SPeter Brune ynorm = &linesearch->ynorm; 66bbd5d0b3SPeter Brune 679105210eSPeter Brune if (!snes->jacobian) { 689105210eSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 699105210eSPeter Brune } 709105210eSPeter Brune 71ea630c6eSPeter Brune /* 72ea630c6eSPeter Brune 73d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 74ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 75ea630c6eSPeter Brune */ 76d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);CHKERRQ(ierr); 77ea630c6eSPeter Brune ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 78ea630c6eSPeter Brune ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 79ea630c6eSPeter Brune ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 80ea630c6eSPeter Brune alpha = alpha / ptAp; 819105210eSPeter Brune ierr = VecAXPY(X,-alpha,Y);CHKERRQ(ierr); 82bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 83bbd5d0b3SPeter Brune 84bbd5d0b3SPeter Brune ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 85bbd5d0b3SPeter Brune ierr = VecNorm(X,NORM_2,xnorm);CHKERRQ(ierr); 86bbd5d0b3SPeter Brune ierr = VecNorm(Y,NORM_2,ynorm);CHKERRQ(ierr); 87ea630c6eSPeter Brune PetscFunctionReturn(0); 88ea630c6eSPeter Brune } 89ea630c6eSPeter Brune 90*b5badacbSBarry Smith /*MC 91*b5badacbSBarry Smith SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG 92*b5badacbSBarry Smith 93*b5badacbSBarry Smith 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*b5badacbSBarry Smith 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*b5badacbSBarry Smith 96*b5badacbSBarry Smith Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian 97*b5badacbSBarry Smith 98*b5badacbSBarry Smith This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone. 99*b5badacbSBarry Smith 100*b5badacbSBarry Smith Level: advanced 101*b5badacbSBarry Smith 102*b5badacbSBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 103*b5badacbSBarry Smith M*/ 104*b5badacbSBarry Smith 1058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 106bbd5d0b3SPeter Brune { 107bbd5d0b3SPeter Brune PetscFunctionBegin; 108f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 1090298fd71SBarry Smith linesearch->ops->destroy = NULL; 1100298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1110298fd71SBarry Smith linesearch->ops->reset = NULL; 1120298fd71SBarry Smith linesearch->ops->view = NULL; 1130298fd71SBarry Smith linesearch->ops->setup = NULL; 114bbd5d0b3SPeter Brune PetscFunctionReturn(0); 115bbd5d0b3SPeter Brune } 116bbd5d0b3SPeter Brune 11788d374b2SPeter Brune /* 118*b5badacbSBarry Smith SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 119*b5badacbSBarry Smith 120*b5badacbSBarry Smith Input Parameter: 121*b5badacbSBarry Smith . snes - the SNES context 122*b5badacbSBarry Smith 123*b5badacbSBarry Smith Application Interface Routine: SNESSetFromOptions() 124*b5badacbSBarry Smith */ 125*b5badacbSBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes) 126*b5badacbSBarry Smith { 127*b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG*)snes->data; 128*b5badacbSBarry Smith PetscErrorCode ierr; 129*b5badacbSBarry Smith PetscBool debug = PETSC_FALSE; 130*b5badacbSBarry Smith SNESNCGType ncgtype=ncg->type; 131*b5badacbSBarry Smith SNESLineSearch linesearch; 132*b5badacbSBarry Smith 133*b5badacbSBarry Smith PetscFunctionBegin; 134*b5badacbSBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES NCG options");CHKERRQ(ierr); 135*b5badacbSBarry Smith 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*b5badacbSBarry Smith if (debug) { 137*b5badacbSBarry Smith ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 138*b5badacbSBarry Smith } 139*b5badacbSBarry Smith ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr); 140*b5badacbSBarry Smith ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 141*b5badacbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 142*b5badacbSBarry Smith if (!snes->linesearch) { 143*b5badacbSBarry Smith ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 144*b5badacbSBarry Smith if (!snes->npc) { 145*b5badacbSBarry Smith ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 146*b5badacbSBarry Smith } else { 147*b5badacbSBarry Smith ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 148*b5badacbSBarry Smith } 149*b5badacbSBarry Smith } 150*b5badacbSBarry Smith PetscFunctionReturn(0); 151*b5badacbSBarry Smith } 152*b5badacbSBarry Smith 153*b5badacbSBarry Smith /* 154*b5badacbSBarry Smith SNESView_NCG - Prints info from the SNESNCG data structure. 155*b5badacbSBarry Smith 156*b5badacbSBarry Smith Input Parameters: 157*b5badacbSBarry Smith + SNES - the SNES context 158*b5badacbSBarry Smith - viewer - visualization context 159*b5badacbSBarry Smith 160*b5badacbSBarry Smith Application Interface Routine: SNESView() 161*b5badacbSBarry Smith */ 162*b5badacbSBarry Smith static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 163*b5badacbSBarry Smith { 164*b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *) snes->data; 165*b5badacbSBarry Smith PetscBool iascii; 166*b5badacbSBarry Smith PetscErrorCode ierr; 167*b5badacbSBarry Smith 168*b5badacbSBarry Smith PetscFunctionBegin; 169*b5badacbSBarry Smith ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 170*b5badacbSBarry Smith if (iascii) { 171*b5badacbSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type]);CHKERRQ(ierr); 172*b5badacbSBarry Smith } 173*b5badacbSBarry Smith PetscFunctionReturn(0); 174*b5badacbSBarry Smith } 175*b5badacbSBarry Smith 176*b5badacbSBarry Smith /* 17788d374b2SPeter Brune 17888d374b2SPeter Brune Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 17988d374b2SPeter Brune 180*b5badacbSBarry Smith This routine is not currently used. I don't know what its intended purpose is. 181*b5badacbSBarry Smith 182*b5badacbSBarry Smith Note it has a hardwired differencing parameter of 1e-5 183*b5badacbSBarry Smith 18488d374b2SPeter Brune */ 18567392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) 18667392de3SBarry Smith { 18788d374b2SPeter Brune PetscErrorCode ierr; 18888d374b2SPeter Brune PetscScalar ftf, ftg, fty, h; 18967392de3SBarry Smith 19088d374b2SPeter Brune PetscFunctionBegin; 19188d374b2SPeter Brune ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 19288d374b2SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 19388d374b2SPeter Brune h = 1e-5*fty / fty; 19488d374b2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 19588d374b2SPeter Brune ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 19688d374b2SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 19788d374b2SPeter Brune ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 19888d374b2SPeter Brune *ytJtf = (ftg - ftf) / h; 19988d374b2SPeter Brune PetscFunctionReturn(0); 20088d374b2SPeter Brune } 20188d374b2SPeter Brune 2020a844d1aSPeter Brune /*@ 2030a844d1aSPeter Brune SNESNCGSetType - Sets the conjugate update type for SNESNCG. 2040a844d1aSPeter Brune 2050a844d1aSPeter Brune Logically Collective on SNES 2060a844d1aSPeter Brune 2070a844d1aSPeter Brune Input Parameters: 2080a844d1aSPeter Brune + snes - the iterative context 2090a844d1aSPeter Brune - btype - update type 2100a844d1aSPeter Brune 2110a844d1aSPeter Brune Options Database: 212*b5badacbSBarry Smith . -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta 2130a844d1aSPeter Brune 2140a844d1aSPeter Brune Level: intermediate 2150a844d1aSPeter Brune 2160a844d1aSPeter Brune SNESNCGSelectTypes: 2170a844d1aSPeter Brune + SNES_NCG_FR - Fletcher-Reeves update 2180a844d1aSPeter Brune . SNES_NCG_PRP - Polak-Ribiere-Polyak update 2190a844d1aSPeter Brune . SNES_NCG_HS - Hestenes-Steifel update 2200a844d1aSPeter Brune . SNES_NCG_DY - Dai-Yuan update 2210a844d1aSPeter Brune - SNES_NCG_CD - Conjugate Descent update 2220a844d1aSPeter Brune 2230a844d1aSPeter Brune Notes: 224*b5badacbSBarry Smith SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions. 225*b5badacbSBarry Smith 226*b5badacbSBarry Smith It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner, 227*b5badacbSBarry Smith that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()? 2280a844d1aSPeter Brune 2290a844d1aSPeter Brune @*/ 2300adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 2310adebc6cSBarry Smith { 2320a844d1aSPeter Brune PetscErrorCode ierr; 2336e111a19SKarl Rupp 2340a844d1aSPeter Brune PetscFunctionBegin; 2350a844d1aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2360a844d1aSPeter Brune ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr); 2370a844d1aSPeter Brune PetscFunctionReturn(0); 2380a844d1aSPeter Brune } 2390a844d1aSPeter Brune 240*b5badacbSBarry Smith static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 2410adebc6cSBarry Smith { 2420a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 2436e111a19SKarl Rupp 2440a844d1aSPeter Brune PetscFunctionBegin; 2450a844d1aSPeter Brune ncg->type = btype; 2460a844d1aSPeter Brune PetscFunctionReturn(0); 2470a844d1aSPeter Brune } 2480a844d1aSPeter Brune 249fef7b6d8SPeter Brune /* 250fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 251fef7b6d8SPeter Brune 252fef7b6d8SPeter Brune Input Parameters: 253fef7b6d8SPeter Brune . snes - the SNES context 254fef7b6d8SPeter Brune 255fef7b6d8SPeter Brune Output Parameter: 256fef7b6d8SPeter Brune . outits - number of iterations until termination 257fef7b6d8SPeter Brune 258fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 259fef7b6d8SPeter Brune */ 260*b5badacbSBarry Smith static PetscErrorCode SNESSolve_NCG(SNES snes) 261fef7b6d8SPeter Brune { 262dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 26332cc618eSPeter Brune Vec X,dX,lX,F,dXold; 264bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0; 26532cc618eSPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 266fef7b6d8SPeter Brune PetscInt maxits, i; 267fef7b6d8SPeter Brune PetscErrorCode ierr; 268422a814eSBarry Smith SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 269f1c6b773SPeter Brune SNESLineSearch linesearch; 270b7281c8aSPeter Brune SNESConvergedReason reason; 27188d374b2SPeter Brune 272fef7b6d8SPeter Brune PetscFunctionBegin; 2736c4ed002SBarry 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); 274c579b300SPatrick Farrell 275fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 276fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 277fef7b6d8SPeter Brune 278fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 279fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 28032cc618eSPeter Brune dXold = snes->work[0]; /* The previous iterate of X */ 28188d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */ 282169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 283fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 284fef7b6d8SPeter Brune 2857601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 2869e764e56SPeter Brune 287e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 288fef7b6d8SPeter Brune snes->iter = 0; 289fef7b6d8SPeter Brune snes->norm = 0.; 290e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 291fef7b6d8SPeter Brune 292bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */ 293a71552e2SPeter Brune 294efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 295be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 296efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 297b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 298b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 299b7281c8aSPeter Brune PetscFunctionReturn(0); 300b7281c8aSPeter Brune } 301a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 302a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 303a71552e2SPeter Brune } else { 304e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 305fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3061aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 307c1c75074SPeter Brune 308fef7b6d8SPeter Brune /* convergence test */ 309a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 310422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 311a71552e2SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 312a71552e2SPeter Brune } 313efd4aadfSBarry Smith if (snes->npc) { 314a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 315be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 316efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 317b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 318b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 319b7281c8aSPeter Brune PetscFunctionReturn(0); 320b7281c8aSPeter Brune } 321a71552e2SPeter Brune } 322a71552e2SPeter Brune } 323a71552e2SPeter Brune ierr = VecCopy(dX,lX);CHKERRQ(ierr); 32432cc618eSPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 325a71552e2SPeter Brune 326a71552e2SPeter Brune 327e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 328fef7b6d8SPeter Brune snes->norm = fnorm; 329e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 330a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 331fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 332fef7b6d8SPeter Brune 333fef7b6d8SPeter Brune /* test convergence */ 334fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 335fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 336fef7b6d8SPeter Brune 337fef7b6d8SPeter Brune /* Call general purpose update function */ 338fef7b6d8SPeter Brune if (snes->ops->update) { 339fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 340fef7b6d8SPeter Brune } 341fef7b6d8SPeter Brune 342fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 343bbd5d0b3SPeter Brune 34409c08436SPeter Brune for (i = 1; i < maxits + 1; i++) { 34529c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 3460a844d1aSPeter Brune if (ncg->type != SNES_NCG_FR) { 34732cc618eSPeter Brune ierr = VecCopy(dX, dXold);CHKERRQ(ierr); 34888d374b2SPeter Brune } 349f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr); 350422a814eSBarry Smith ierr = SNESLineSearchGetReason(linesearch, &lsresult);CHKERRQ(ierr); 351422a814eSBarry Smith ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 352422a814eSBarry Smith if (lsresult) { 353fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 354fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 3555d10c001SPeter Brune PetscFunctionReturn(0); 356fef7b6d8SPeter Brune } 357fef7b6d8SPeter Brune } 358e71169deSBarry Smith if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 359fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 3605d10c001SPeter Brune PetscFunctionReturn(0); 361fef7b6d8SPeter Brune } 362fef7b6d8SPeter Brune /* Monitor convergence */ 363e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 364169e2be8SPeter Brune snes->iter = i; 365fef7b6d8SPeter Brune snes->norm = fnorm; 366c1e67a49SFande Kong snes->xnorm = xnorm; 367c1e67a49SFande Kong snes->ynorm = ynorm; 368e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 369a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 370fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 371302dbbaeSPeter Brune 372fef7b6d8SPeter Brune /* Test for convergence */ 373d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3745d10c001SPeter Brune if (snes->reason) PetscFunctionReturn(0); 375302dbbaeSPeter Brune 376302dbbaeSPeter Brune /* Call general purpose update function */ 377302dbbaeSPeter Brune if (snes->ops->update) { 378302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 379302dbbaeSPeter Brune } 380efd4aadfSBarry Smith if (snes->npc) { 381a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 382be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr); 383efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 384b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 385b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 386b7281c8aSPeter Brune PetscFunctionReturn(0); 387b7281c8aSPeter Brune } 388a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 389a71552e2SPeter Brune } else { 390be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr); 391efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 392b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 393b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 394b7281c8aSPeter Brune PetscFunctionReturn(0); 395b7281c8aSPeter Brune } 396302dbbaeSPeter Brune } 397a3685310SPeter Brune } else { 398a3685310SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 399302dbbaeSPeter Brune } 400302dbbaeSPeter Brune 40129c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 4020a844d1aSPeter Brune switch (ncg->type) { 4030a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */ 40432cc618eSPeter Brune dXolddotdXold= dXdotdX; 40532cc618eSPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 40632cc618eSPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold); 407dfb256c7SPeter Brune break; 4080a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 40932cc618eSPeter Brune dXolddotdXold= dXdotdX; 41032cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr); 41132cc618eSPeter Brune ierr = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 41232cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 41332cc618eSPeter Brune ierr = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 41432cc618eSPeter Brune beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 415dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 416dfb256c7SPeter Brune break; 4170a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */ 41832cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 41932cc618eSPeter Brune ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 42032cc618eSPeter Brune ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 42132cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 42232cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 42332cc618eSPeter Brune ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 42432cc618eSPeter Brune ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 42532cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 42632cc618eSPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 427dfb256c7SPeter Brune break; 4280a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */ 42932cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 43032cc618eSPeter Brune ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 43132cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 43232cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 43332cc618eSPeter Brune ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 43432cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 43532cc618eSPeter Brune beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr); 436dfb256c7SPeter Brune break; 4370a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */ 43832cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 43932cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 44032cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 44132cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 44232cc618eSPeter Brune beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr); 443dfb256c7SPeter Brune break; 444dfb256c7SPeter Brune } 445dfb256c7SPeter Brune if (ncg->monitor) { 4467904a332SBarry Smith ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);CHKERRQ(ierr); 447dfb256c7SPeter Brune } 448dfb256c7SPeter Brune ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 449fef7b6d8SPeter Brune } 450fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 451fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 452fef7b6d8SPeter Brune PetscFunctionReturn(0); 453fef7b6d8SPeter Brune } 454fef7b6d8SPeter Brune 455fef7b6d8SPeter Brune /*MC 456fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 457fef7b6d8SPeter Brune 458fef7b6d8SPeter Brune Level: beginner 459fef7b6d8SPeter Brune 460fef7b6d8SPeter Brune Options Database: 4614f02bc6aSBarry Smith + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp. 46241a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type. 463*b5badacbSBarry Smith - -snes_ncg_monitor - Print the beta values used in the ncg iteration, . 464fef7b6d8SPeter Brune 46595452b02SPatrick Sanan Notes: 46695452b02SPatrick Sanan This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 467fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 468fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 469fef7b6d8SPeter Brune 4702d547940SBarry Smith Only supports left non-linear preconditioning. 4712d547940SBarry Smith 472*b5badacbSBarry Smith Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then 473*b5badacbSBarry Smith SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR 474*b5badacbSBarry Smith 4754f02bc6aSBarry Smith References: 47696a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers", 4774f02bc6aSBarry Smith SIAM Review, 57(4), 2015 478cc3c3c0fSMatthew G. Knepley 479fef7b6d8SPeter Brune 480*b5badacbSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType() 481fef7b6d8SPeter Brune M*/ 4828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 483fef7b6d8SPeter Brune { 484fef7b6d8SPeter Brune PetscErrorCode ierr; 485ea630c6eSPeter Brune SNES_NCG * neP; 486fef7b6d8SPeter Brune 487fef7b6d8SPeter Brune PetscFunctionBegin; 488fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 489fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 490fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 491fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 492fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 493fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 494fef7b6d8SPeter Brune 495fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 496efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 497efd4aadfSBarry Smith snes->npcside = PC_LEFT; 498fef7b6d8SPeter Brune 4994fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 5004fc747eaSLawrence Mitchell 50188976e71SPeter Brune if (!snes->tolerancesset) { 5020e444f03SPeter Brune snes->max_funcs = 30000; 5030e444f03SPeter Brune snes->max_its = 10000; 504c522fa08SPeter Brune snes->stol = 1e-20; 50588976e71SPeter Brune } 5060e444f03SPeter Brune 507b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 508fef7b6d8SPeter Brune snes->data = (void*) neP; 5090298fd71SBarry Smith neP->monitor = NULL; 5100a844d1aSPeter Brune neP->type = SNES_NCG_PRP; 511bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr); 512fef7b6d8SPeter Brune PetscFunctionReturn(0); 513fef7b6d8SPeter Brune } 514