10a844d1aSPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/ 29e5d0892SLisandro Dalcin const char *const SNESNCGTypes[] = {"FR", "PRP", "HS", "DY", "CD", "SNESNCGType", "SNES_NCG_", NULL}; 3fef7b6d8SPeter Brune 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NCG(SNES snes) 5d71ae5a4SJacob Faibussowitsch { 6fef7b6d8SPeter Brune PetscFunctionBegin; 72e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL)); 89566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 93ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10fef7b6d8SPeter Brune } 11fef7b6d8SPeter Brune 12d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NCG(SNES snes) 13d71ae5a4SJacob Faibussowitsch { 14fef7b6d8SPeter Brune PetscFunctionBegin; 159566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 2)); 1608401ef6SPierre Jolivet PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 176c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19fef7b6d8SPeter Brune } 20fef7b6d8SPeter Brune 21d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 22d71ae5a4SJacob Faibussowitsch { 23ea630c6eSPeter Brune PetscScalar alpha, ptAp; 24bbd5d0b3SPeter Brune Vec X, Y, F, W; 25bbd5d0b3SPeter Brune SNES snes; 26bbd5d0b3SPeter Brune PetscReal *fnorm, *xnorm, *ynorm; 27ea630c6eSPeter Brune 28ea630c6eSPeter Brune PetscFunctionBegin; 299566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 30bbd5d0b3SPeter Brune X = linesearch->vec_sol; 31bbd5d0b3SPeter Brune W = linesearch->vec_sol_new; 32bbd5d0b3SPeter Brune F = linesearch->vec_func; 33bbd5d0b3SPeter Brune Y = linesearch->vec_update; 34bbd5d0b3SPeter Brune fnorm = &linesearch->fnorm; 35bbd5d0b3SPeter Brune xnorm = &linesearch->xnorm; 36bbd5d0b3SPeter Brune ynorm = &linesearch->ynorm; 37bbd5d0b3SPeter Brune 3848a46eb9SPierre Jolivet if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes)); 399105210eSPeter Brune 40ea630c6eSPeter Brune /* 41d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 42ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 43ea630c6eSPeter Brune */ 449566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 459566063dSJacob Faibussowitsch PetscCall(VecDot(F, F, &alpha)); 469566063dSJacob Faibussowitsch PetscCall(MatMult(snes->jacobian, Y, W)); 479566063dSJacob Faibussowitsch PetscCall(VecDot(Y, W, &ptAp)); 48ea630c6eSPeter Brune alpha = alpha / ptAp; 499566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, -alpha, Y)); 509566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 51bbd5d0b3SPeter Brune 529566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 539566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, xnorm)); 549566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, ynorm)); 553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56ea630c6eSPeter Brune } 57ea630c6eSPeter Brune 58b5badacbSBarry Smith /*MC 59f6dfbefdSBarry Smith SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG` 60b5badacbSBarry Smith 61b5badacbSBarry 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. 62b5badacbSBarry 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. 63b5badacbSBarry Smith 64420bcc1bSBarry Smith Level: advanced 65420bcc1bSBarry Smith 66f6dfbefdSBarry Smith Notes: 67f6dfbefdSBarry Smith This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian 68b5badacbSBarry Smith 69b5badacbSBarry Smith This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone. 70b5badacbSBarry Smith 710241b274SPierre Jolivet .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()` 72b5badacbSBarry Smith M*/ 73b5badacbSBarry Smith 74d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 75d71ae5a4SJacob Faibussowitsch { 76bbd5d0b3SPeter Brune PetscFunctionBegin; 77f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 780298fd71SBarry Smith linesearch->ops->destroy = NULL; 790298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 800298fd71SBarry Smith linesearch->ops->reset = NULL; 810298fd71SBarry Smith linesearch->ops->view = NULL; 820298fd71SBarry Smith linesearch->ops->setup = NULL; 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84bbd5d0b3SPeter Brune } 85bbd5d0b3SPeter Brune 86ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems PetscOptionsObject) 87d71ae5a4SJacob Faibussowitsch { 88b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *)snes->data; 89b5badacbSBarry Smith PetscBool debug = PETSC_FALSE; 90b5badacbSBarry Smith SNESNCGType ncgtype = ncg->type; 91b5badacbSBarry Smith SNESLineSearch linesearch; 92b5badacbSBarry Smith 93b5badacbSBarry Smith PetscFunctionBegin; 94d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options"); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL)); 96ad540459SPierre Jolivet if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 979566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL)); 989566063dSJacob Faibussowitsch PetscCall(SNESNCGSetType(snes, ncgtype)); 99d0609cedSBarry Smith PetscOptionsHeadEnd(); 100b5badacbSBarry Smith if (!snes->linesearch) { 1019566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 102ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 103b5badacbSBarry Smith if (!snes->npc) { 1049566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP)); 105b5badacbSBarry Smith } else { 106a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT)); 107b5badacbSBarry Smith } 108b5badacbSBarry Smith } 109ec786807SJed Brown } 1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 111b5badacbSBarry Smith } 112b5badacbSBarry Smith 113d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 114d71ae5a4SJacob Faibussowitsch { 115b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *)snes->data; 116*9f196a02SMartin Diehl PetscBool isascii; 117b5badacbSBarry Smith 118b5badacbSBarry Smith PetscFunctionBegin; 119*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 120*9f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type])); 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122b5badacbSBarry Smith } 123b5badacbSBarry Smith 1240a844d1aSPeter Brune /*@ 125f6dfbefdSBarry Smith SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`. 1260a844d1aSPeter Brune 127c3339decSBarry Smith Logically Collective 1280a844d1aSPeter Brune 1290a844d1aSPeter Brune Input Parameters: 1300a844d1aSPeter Brune + snes - the iterative context 131ceaaa498SBarry Smith - btype - update type, see `SNESNCGType` 1320a844d1aSPeter Brune 1333c7db156SBarry Smith Options Database Key: 134b5badacbSBarry Smith . -snes_ncg_type <prp,fr,hs,dy,cd> - strategy for selecting algorithm for computing beta 1350a844d1aSPeter Brune 1360a844d1aSPeter Brune Level: intermediate 1370a844d1aSPeter Brune 1380a844d1aSPeter Brune Notes: 139f6dfbefdSBarry Smith `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions. 140b5badacbSBarry Smith 141b5badacbSBarry Smith It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner, 142f6dfbefdSBarry Smith that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`? 1430a844d1aSPeter Brune 144420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD` 1450a844d1aSPeter Brune @*/ 146d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 147d71ae5a4SJacob Faibussowitsch { 1480a844d1aSPeter Brune PetscFunctionBegin; 1490a844d1aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 150cac4c232SBarry Smith PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype)); 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1520a844d1aSPeter Brune } 1530a844d1aSPeter Brune 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 155d71ae5a4SJacob Faibussowitsch { 1560a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data; 1576e111a19SKarl Rupp 1580a844d1aSPeter Brune PetscFunctionBegin; 1590a844d1aSPeter Brune ncg->type = btype; 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1610a844d1aSPeter Brune } 1620a844d1aSPeter Brune 163fef7b6d8SPeter Brune /* 164fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 165fef7b6d8SPeter Brune 1662fe279fdSBarry Smith Input Parameter: 1672fe279fdSBarry Smith . snes - the `SNES` context 168fef7b6d8SPeter Brune 169fef7b6d8SPeter Brune Output Parameter: 170fef7b6d8SPeter Brune . outits - number of iterations until termination 171fef7b6d8SPeter Brune 172fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 173fef7b6d8SPeter Brune */ 174d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NCG(SNES snes) 175d71ae5a4SJacob Faibussowitsch { 176dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data; 17732cc618eSPeter Brune Vec X, dX, lX, F, dXold; 178bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0; 17932cc618eSPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 180fef7b6d8SPeter Brune PetscInt maxits, i; 181422a814eSBarry Smith SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 182f1c6b773SPeter Brune SNESLineSearch linesearch; 183b7281c8aSPeter Brune SNESConvergedReason reason; 18488d374b2SPeter Brune 185fef7b6d8SPeter Brune PetscFunctionBegin; 1860b121fc5SBarry Smith PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 187c579b300SPatrick Farrell 1889566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 189fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 190fef7b6d8SPeter Brune 191fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 192fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 19332cc618eSPeter Brune dXold = snes->work[0]; /* The previous iterate of X */ 19488d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */ 195169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 196fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 197fef7b6d8SPeter Brune 1989566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 1999e764e56SPeter Brune 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 201fef7b6d8SPeter Brune snes->iter = 0; 202fef7b6d8SPeter Brune snes->norm = 0.; 2039566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 204fef7b6d8SPeter Brune 205bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */ 206a71552e2SPeter Brune 207efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 2089566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, dX)); 2099566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 210b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 211b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 213b7281c8aSPeter Brune } 2149566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, F)); 2159566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 216a71552e2SPeter Brune } else { 217e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2189566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 2191aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 220c1c75074SPeter Brune 221fef7b6d8SPeter Brune /* convergence test */ 2229566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 223422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 2249566063dSJacob Faibussowitsch PetscCall(VecCopy(F, dX)); 225a71552e2SPeter Brune } 226efd4aadfSBarry Smith if (snes->npc) { 227a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2289566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, dX)); 2299566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 230b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 231b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233b7281c8aSPeter Brune } 234a71552e2SPeter Brune } 235a71552e2SPeter Brune } 2369566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, lX)); 2379566063dSJacob Faibussowitsch PetscCall(VecDot(dX, dX, &dXdotdX)); 238a71552e2SPeter Brune 2399566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 240fef7b6d8SPeter Brune snes->norm = fnorm; 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2429566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 243fef7b6d8SPeter Brune 244fef7b6d8SPeter Brune /* test convergence */ 2452d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 2462d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 2473ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 248fef7b6d8SPeter Brune 249fef7b6d8SPeter Brune /* Call general purpose update function */ 250dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 251fef7b6d8SPeter Brune 252fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 253bbd5d0b3SPeter Brune 25409c08436SPeter Brune for (i = 1; i < maxits + 1; i++) { 25529c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 25648a46eb9SPierre Jolivet if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold)); 2579566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX)); 2589566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lsresult)); 2599566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 260422a814eSBarry Smith if (lsresult) { 261fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 262fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 264fef7b6d8SPeter Brune } 265fef7b6d8SPeter Brune } 266e71169deSBarry Smith if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 267fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269fef7b6d8SPeter Brune } 270fef7b6d8SPeter Brune /* Monitor convergence */ 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 272169e2be8SPeter Brune snes->iter = i; 273fef7b6d8SPeter Brune snes->norm = fnorm; 274c1e67a49SFande Kong snes->xnorm = xnorm; 275c1e67a49SFande Kong snes->ynorm = ynorm; 2769566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2779566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 278302dbbaeSPeter Brune 279fef7b6d8SPeter Brune /* Test for convergence */ 2802d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 2812d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2823ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 283302dbbaeSPeter Brune 284302dbbaeSPeter Brune /* Call general purpose update function */ 285dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 286efd4aadfSBarry Smith if (snes->npc) { 287a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 2889566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, dX)); 2899566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 290b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 291b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 293b7281c8aSPeter Brune } 2949566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, F)); 295a71552e2SPeter Brune } else { 2969566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, dX)); 2979566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 298b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 299b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 301b7281c8aSPeter Brune } 302302dbbaeSPeter Brune } 303a3685310SPeter Brune } else { 3049566063dSJacob Faibussowitsch PetscCall(VecCopy(F, dX)); 305302dbbaeSPeter Brune } 306302dbbaeSPeter Brune 30729c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 3080a844d1aSPeter Brune switch (ncg->type) { 3090a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */ 31032cc618eSPeter Brune dXolddotdXold = dXdotdX; 3119566063dSJacob Faibussowitsch PetscCall(VecDot(dX, dX, &dXdotdX)); 31232cc618eSPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold); 313dfb256c7SPeter Brune break; 3140a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 31532cc618eSPeter Brune dXolddotdXold = dXdotdX; 316d3981886SPierre Jolivet PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3179566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dXold, dX, &dXdotdXold)); 3189566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3199566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dXold, dX, &dXdotdXold)); 320f4f49eeaSPierre Jolivet beta = PetscRealPart((dXdotdX - dXdotdXold) / dXolddotdXold); 321dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 322dfb256c7SPeter Brune break; 3230a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */ 3249566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3259566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dXold, &dXdotdXold)); 3269566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dX, &lXdotdX)); 3279566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 3289566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3299566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dXold, &dXdotdXold)); 3309566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dX, &lXdotdX)); 3319566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 33232cc618eSPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 333dfb256c7SPeter Brune break; 3340a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */ 3359566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3369566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dX, &lXdotdX)); 3379566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 3389566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3399566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dX, &lXdotdX)); 3409566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 3412f613bf5SBarry Smith beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX)); 342dfb256c7SPeter Brune break; 3430a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */ 3449566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3459566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 3469566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3479566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 3482f613bf5SBarry Smith beta = PetscRealPart(dXdotdX / lXdotdXold); 349dfb256c7SPeter Brune break; 350dfb256c7SPeter Brune } 35148a46eb9SPierre Jolivet if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta)); 3529566063dSJacob Faibussowitsch PetscCall(VecAYPX(lX, beta, dX)); 353fef7b6d8SPeter Brune } 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355fef7b6d8SPeter Brune } 356fef7b6d8SPeter Brune 357fef7b6d8SPeter Brune /*MC 3581d27aa22SBarry Smith SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems {cite}`bruneknepleysmithtu15`. 359fef7b6d8SPeter Brune 360fef7b6d8SPeter Brune Level: beginner 361fef7b6d8SPeter Brune 362f6dfbefdSBarry Smith Options Database Keys: 36362842358SBarry Smith + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is `prp`. 36441a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type. 3651d27aa22SBarry Smith - -snes_ncg_monitor - Print the beta values nonlinear Conjugate-Gradient used in the iteration, . 366fef7b6d8SPeter Brune 36795452b02SPatrick Sanan Notes: 3681d27aa22SBarry Smith This solves the nonlinear system of equations $ F(x) = 0 $ using the nonlinear generalization of the conjugate 369fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 3701d27aa22SBarry Smith chooses the initial search direction as $ F(x) $ for the initial guess $x$. 371fef7b6d8SPeter Brune 3722d547940SBarry Smith Only supports left non-linear preconditioning. 3732d547940SBarry Smith 37462842358SBarry Smith Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with `-npc_snes_type` <type>, `SNESSetNPC()`, or `SNESGetNPC()` then 375a99ef635SJonas Heinzmann `SNESLINESEARCHSECANT` is used. Also supports the special-purpose line search `SNESLINESEARCHNCGLINEAR` 376b5badacbSBarry Smith 377420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()` 378fef7b6d8SPeter Brune M*/ 379d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 380d71ae5a4SJacob Faibussowitsch { 381ea630c6eSPeter Brune SNES_NCG *neP; 382fef7b6d8SPeter Brune 383fef7b6d8SPeter Brune PetscFunctionBegin; 384fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 385fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 386fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 387fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 388fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 389fef7b6d8SPeter Brune 390fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 391efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 392efd4aadfSBarry Smith snes->npcside = PC_LEFT; 393fef7b6d8SPeter Brune 3944fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3954fc747eaSLawrence Mitchell 39677e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 39777e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_funcs, 30000); 39877e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_its, 10000); 39977e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, stol, 1e-20); 4000e444f03SPeter Brune 4014dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 402fef7b6d8SPeter Brune snes->data = (void *)neP; 4030298fd71SBarry Smith neP->monitor = NULL; 4040a844d1aSPeter Brune neP->type = SNES_NCG_PRP; 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG)); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407fef7b6d8SPeter Brune } 408