xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 2fe279fdf3e687a416e4eadb7d3c7a82d60442c6)
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 SNESReset_NCG(SNES snes)
5d71ae5a4SJacob Faibussowitsch {
6fef7b6d8SPeter Brune   PetscFunctionBegin;
73ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8fef7b6d8SPeter Brune }
9fef7b6d8SPeter Brune 
10d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NCG(SNES snes)
11d71ae5a4SJacob Faibussowitsch {
12fef7b6d8SPeter Brune   PetscFunctionBegin;
132e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL));
149566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16fef7b6d8SPeter Brune }
17fef7b6d8SPeter Brune 
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NCG(SNES snes)
19d71ae5a4SJacob Faibussowitsch {
20fef7b6d8SPeter Brune   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 2));
2208401ef6SPierre Jolivet   PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
236c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25fef7b6d8SPeter Brune }
26fef7b6d8SPeter Brune 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
28d71ae5a4SJacob Faibussowitsch {
29ea630c6eSPeter Brune   PetscScalar alpha, ptAp;
30bbd5d0b3SPeter Brune   Vec         X, Y, F, W;
31bbd5d0b3SPeter Brune   SNES        snes;
32bbd5d0b3SPeter Brune   PetscReal  *fnorm, *xnorm, *ynorm;
33ea630c6eSPeter Brune 
34ea630c6eSPeter Brune   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
36bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
37bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
38bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
39bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
40bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
41bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
42bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
43bbd5d0b3SPeter Brune 
4448a46eb9SPierre Jolivet   if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));
459105210eSPeter Brune 
46ea630c6eSPeter Brune   /*
47d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
48ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
49ea630c6eSPeter Brune    */
509566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
519566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &alpha));
529566063dSJacob Faibussowitsch   PetscCall(MatMult(snes->jacobian, Y, W));
539566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, W, &ptAp));
54ea630c6eSPeter Brune   alpha = alpha / ptAp;
559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, -alpha, Y));
569566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, F));
57bbd5d0b3SPeter Brune 
589566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, fnorm));
599566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, xnorm));
609566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, ynorm));
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62ea630c6eSPeter Brune }
63ea630c6eSPeter Brune 
64b5badacbSBarry Smith /*MC
65f6dfbefdSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG`
66b5badacbSBarry Smith 
67b5badacbSBarry 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.
68b5badacbSBarry 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.
69b5badacbSBarry Smith 
70f6dfbefdSBarry Smith    Notes:
71f6dfbefdSBarry Smith    This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
72b5badacbSBarry Smith 
73b5badacbSBarry Smith    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
74b5badacbSBarry Smith 
75b5badacbSBarry Smith    Level: advanced
76b5badacbSBarry Smith 
77db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
78b5badacbSBarry Smith M*/
79b5badacbSBarry Smith 
80d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
81d71ae5a4SJacob Faibussowitsch {
82bbd5d0b3SPeter Brune   PetscFunctionBegin;
83f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
840298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
850298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
860298fd71SBarry Smith   linesearch->ops->reset          = NULL;
870298fd71SBarry Smith   linesearch->ops->view           = NULL;
880298fd71SBarry Smith   linesearch->ops->setup          = NULL;
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90bbd5d0b3SPeter Brune }
91bbd5d0b3SPeter Brune 
92d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject)
93d71ae5a4SJacob Faibussowitsch {
94b5badacbSBarry Smith   SNES_NCG      *ncg     = (SNES_NCG *)snes->data;
95b5badacbSBarry Smith   PetscBool      debug   = PETSC_FALSE;
96b5badacbSBarry Smith   SNESNCGType    ncgtype = ncg->type;
97b5badacbSBarry Smith   SNESLineSearch linesearch;
98b5badacbSBarry Smith 
99b5badacbSBarry Smith   PetscFunctionBegin;
100d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
1019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
102ad540459SPierre Jolivet   if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(SNESNCGSetType(snes, ncgtype));
105d0609cedSBarry Smith   PetscOptionsHeadEnd();
106b5badacbSBarry Smith   if (!snes->linesearch) {
1079566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
108ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
109b5badacbSBarry Smith       if (!snes->npc) {
1109566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
111b5badacbSBarry Smith       } else {
1129566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
113b5badacbSBarry Smith       }
114b5badacbSBarry Smith     }
115ec786807SJed Brown   }
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117b5badacbSBarry Smith }
118b5badacbSBarry Smith 
119d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
120d71ae5a4SJacob Faibussowitsch {
121b5badacbSBarry Smith   SNES_NCG *ncg = (SNES_NCG *)snes->data;
122b5badacbSBarry Smith   PetscBool iascii;
123b5badacbSBarry Smith 
124b5badacbSBarry Smith   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
12648a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
1273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128b5badacbSBarry Smith }
129b5badacbSBarry Smith 
130b5badacbSBarry Smith /*
13188d374b2SPeter Brune 
13288d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
13388d374b2SPeter Brune 
134b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
135b5badacbSBarry Smith 
136b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
137b5badacbSBarry Smith 
13888d374b2SPeter Brune  */
139d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar *ytJtf)
140d71ae5a4SJacob Faibussowitsch {
14188d374b2SPeter Brune   PetscScalar ftf, ftg, fty, h;
14267392de3SBarry Smith 
14388d374b2SPeter Brune   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &ftf));
1459566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty));
14688d374b2SPeter Brune   h = 1e-5 * fty / fty;
1479566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
1489566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W, -h, Y)); /* this is arbitrary */
1499566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, W, G));
1509566063dSJacob Faibussowitsch   PetscCall(VecDot(G, F, &ftg));
15188d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15388d374b2SPeter Brune }
15488d374b2SPeter Brune 
1550a844d1aSPeter Brune /*@
156f6dfbefdSBarry Smith     SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`.
1570a844d1aSPeter Brune 
158c3339decSBarry Smith     Logically Collective
1590a844d1aSPeter Brune 
1600a844d1aSPeter Brune     Input Parameters:
1610a844d1aSPeter Brune +   snes - the iterative context
1620a844d1aSPeter Brune -   btype - update type
1630a844d1aSPeter Brune 
1643c7db156SBarry Smith     Options Database Key:
165b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
1660a844d1aSPeter Brune 
1670a844d1aSPeter Brune     Level: intermediate
1680a844d1aSPeter Brune 
169f6dfbefdSBarry Smith     `SNESNCGType`s:
170f6dfbefdSBarry Smith +   `SNES_NCG_FR` - Fletcher-Reeves update
171f6dfbefdSBarry Smith .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update
172f6dfbefdSBarry Smith .   `SNES_NCG_HS` - Hestenes-Steifel update
173f6dfbefdSBarry Smith .   `SNES_NCG_DY` - Dai-Yuan update
174f6dfbefdSBarry Smith -   `SNES_NCG_CD` - Conjugate Descent update
1750a844d1aSPeter Brune 
1760a844d1aSPeter Brune    Notes:
177f6dfbefdSBarry Smith    `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.
178b5badacbSBarry Smith 
179b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
180f6dfbefdSBarry Smith    that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`?
1810a844d1aSPeter Brune 
182f6dfbefdSBarry Smith    Developer Note:
183f6dfbefdSBarry Smith    There should be a `SNESNCGSetType()`
184f6dfbefdSBarry Smith 
185f6dfbefdSBarry Smith .seealso: `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
1860a844d1aSPeter Brune @*/
187d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
188d71ae5a4SJacob Faibussowitsch {
1890a844d1aSPeter Brune   PetscFunctionBegin;
1900a844d1aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
191cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1930a844d1aSPeter Brune }
1940a844d1aSPeter Brune 
195d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
196d71ae5a4SJacob Faibussowitsch {
1970a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG *)snes->data;
1986e111a19SKarl Rupp 
1990a844d1aSPeter Brune   PetscFunctionBegin;
2000a844d1aSPeter Brune   ncg->type = btype;
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2020a844d1aSPeter Brune }
2030a844d1aSPeter Brune 
204fef7b6d8SPeter Brune /*
205fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
206fef7b6d8SPeter Brune 
207*2fe279fdSBarry Smith   Input Parameter:
208*2fe279fdSBarry Smith . snes - the `SNES` context
209fef7b6d8SPeter Brune 
210fef7b6d8SPeter Brune   Output Parameter:
211fef7b6d8SPeter Brune . outits - number of iterations until termination
212fef7b6d8SPeter Brune 
213fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
214fef7b6d8SPeter Brune */
215d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NCG(SNES snes)
216d71ae5a4SJacob Faibussowitsch {
217dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
21832cc618eSPeter Brune   Vec                  X, dX, lX, F, dXold;
219bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
22032cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
221fef7b6d8SPeter Brune   PetscInt             maxits, i;
222422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
223f1c6b773SPeter Brune   SNESLineSearch       linesearch;
224b7281c8aSPeter Brune   SNESConvergedReason  reason;
22588d374b2SPeter Brune 
226fef7b6d8SPeter Brune   PetscFunctionBegin;
2270b121fc5SBarry 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);
228c579b300SPatrick Farrell 
2299566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
230fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
231fef7b6d8SPeter Brune 
232fef7b6d8SPeter Brune   maxits = snes->max_its;        /* maximum number of iterations */
233fef7b6d8SPeter Brune   X      = snes->vec_sol;        /* X^n */
23432cc618eSPeter Brune   dXold  = snes->work[0];        /* The previous iterate of X */
23588d374b2SPeter Brune   dX     = snes->work[1];        /* the preconditioned direction */
236169e2be8SPeter Brune   lX     = snes->vec_sol_update; /* the conjugate direction */
237fef7b6d8SPeter Brune   F      = snes->vec_func;       /* residual vector */
238fef7b6d8SPeter Brune 
2399566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
2409e764e56SPeter Brune 
2419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
242fef7b6d8SPeter Brune   snes->iter = 0;
243fef7b6d8SPeter Brune   snes->norm = 0.;
2449566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
245fef7b6d8SPeter Brune 
246bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
247a71552e2SPeter Brune 
248efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2499566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2509566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
251b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
252b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
2533ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
254b7281c8aSPeter Brune     }
2559566063dSJacob Faibussowitsch     PetscCall(VecCopy(dX, F));
2569566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
257a71552e2SPeter Brune   } else {
258e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2599566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2601aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
261c1c75074SPeter Brune 
262fef7b6d8SPeter Brune     /* convergence test */
2639566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
264422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2659566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, dX));
266a71552e2SPeter Brune   }
267efd4aadfSBarry Smith   if (snes->npc) {
268a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2699566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, dX));
2709566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
271b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
272b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2733ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
274b7281c8aSPeter Brune       }
275a71552e2SPeter Brune     }
276a71552e2SPeter Brune   }
2779566063dSJacob Faibussowitsch   PetscCall(VecCopy(dX, lX));
2789566063dSJacob Faibussowitsch   PetscCall(VecDot(dX, dX, &dXdotdX));
279a71552e2SPeter Brune 
2809566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
281fef7b6d8SPeter Brune   snes->norm = fnorm;
2829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2839566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
2849566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
285fef7b6d8SPeter Brune 
286fef7b6d8SPeter Brune   /* test convergence */
287dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
2883ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
289fef7b6d8SPeter Brune 
290fef7b6d8SPeter Brune   /* Call general purpose update function */
291dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
292fef7b6d8SPeter Brune 
293fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
294bbd5d0b3SPeter Brune 
29509c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
29629c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
29748a46eb9SPierre Jolivet     if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
2989566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
2999566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
3009566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
301422a814eSBarry Smith     if (lsresult) {
302fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
303fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3043ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
305fef7b6d8SPeter Brune       }
306fef7b6d8SPeter Brune     }
307e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
308fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3093ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
310fef7b6d8SPeter Brune     }
311fef7b6d8SPeter Brune     /* Monitor convergence */
3129566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
313169e2be8SPeter Brune     snes->iter  = i;
314fef7b6d8SPeter Brune     snes->norm  = fnorm;
315c1e67a49SFande Kong     snes->xnorm = xnorm;
316c1e67a49SFande Kong     snes->ynorm = ynorm;
3179566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3189566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
3199566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
320302dbbaeSPeter Brune 
321fef7b6d8SPeter Brune     /* Test for convergence */
322dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
3233ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
324302dbbaeSPeter Brune 
325302dbbaeSPeter Brune     /* Call general purpose update function */
326dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
327efd4aadfSBarry Smith     if (snes->npc) {
328a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3299566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, dX));
3309566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
331b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
332b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
3333ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
334b7281c8aSPeter Brune         }
3359566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX, F));
336a71552e2SPeter Brune       } else {
3379566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, dX));
3389566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
339b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
340b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
3413ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
342b7281c8aSPeter Brune         }
343302dbbaeSPeter Brune       }
344a3685310SPeter Brune     } else {
3459566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, dX));
346302dbbaeSPeter Brune     }
347302dbbaeSPeter Brune 
34829c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3490a844d1aSPeter Brune     switch (ncg->type) {
3500a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
35132cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3529566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
35332cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / dXolddotdXold);
354dfb256c7SPeter Brune       break;
3550a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
35632cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3579566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
3589566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
3599566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3609566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
36132cc618eSPeter Brune       beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
362dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
363dfb256c7SPeter Brune       break;
3640a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
3659566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3669566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
3679566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3689566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3699566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3709566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
3719566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3729566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
37332cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
374dfb256c7SPeter Brune       break;
3750a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
3769566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3779566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3789566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3799566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3809566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3819566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
3822f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
383dfb256c7SPeter Brune       break;
3840a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
3859566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3869566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3879566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3889566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
3892f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
390dfb256c7SPeter Brune       break;
391dfb256c7SPeter Brune     }
39248a46eb9SPierre Jolivet     if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
3939566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
394fef7b6d8SPeter Brune   }
39563a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
396fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398fef7b6d8SPeter Brune }
399fef7b6d8SPeter Brune 
400fef7b6d8SPeter Brune /*MC
401f6dfbefdSBarry Smith   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems.
402fef7b6d8SPeter Brune 
403fef7b6d8SPeter Brune   Level: beginner
404fef7b6d8SPeter Brune 
405f6dfbefdSBarry Smith   Options Database Keys:
4064f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
40741a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
408b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
409fef7b6d8SPeter Brune 
41095452b02SPatrick Sanan    Notes:
41195452b02SPatrick Sanan    This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
412fef7b6d8SPeter Brune    gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
413fef7b6d8SPeter Brune    chooses the initial search direction as F(x) for the initial guess x.
414fef7b6d8SPeter Brune 
4152d547940SBarry Smith    Only supports left non-linear preconditioning.
4162d547940SBarry Smith 
417f6dfbefdSBarry Smith    Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then
418f6dfbefdSBarry Smith    `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR`
419b5badacbSBarry Smith 
4204f02bc6aSBarry Smith    References:
421606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4224f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
423cc3c3c0fSMatthew G. Knepley 
424f6dfbefdSBarry Smith .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
425fef7b6d8SPeter Brune M*/
426d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
427d71ae5a4SJacob Faibussowitsch {
428ea630c6eSPeter Brune   SNES_NCG *neP;
429fef7b6d8SPeter Brune 
430fef7b6d8SPeter Brune   PetscFunctionBegin;
431fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
432fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
433fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
434fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
435fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
436fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
437fef7b6d8SPeter Brune 
438fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
439efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
440efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
441fef7b6d8SPeter Brune 
4424fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4434fc747eaSLawrence Mitchell 
44488976e71SPeter Brune   if (!snes->tolerancesset) {
4450e444f03SPeter Brune     snes->max_funcs = 30000;
4460e444f03SPeter Brune     snes->max_its   = 10000;
447c522fa08SPeter Brune     snes->stol      = 1e-20;
44888976e71SPeter Brune   }
4490e444f03SPeter Brune 
4504dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
451fef7b6d8SPeter Brune   snes->data   = (void *)neP;
4520298fd71SBarry Smith   neP->monitor = NULL;
4530a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
4549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
456fef7b6d8SPeter Brune }
457