xref: /petsc/src/snes/impls/ncg/snesncg.c (revision ceaaa4989964adb3f5eb130cb04b8f49c83e49c9)
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
162*ceaaa498SBarry Smith - btype - update type, see `SNESNCGType`
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 
1690a844d1aSPeter Brune   Notes:
170f6dfbefdSBarry Smith   `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.
171b5badacbSBarry Smith 
172b5badacbSBarry Smith   It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
173f6dfbefdSBarry Smith   that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`?
1740a844d1aSPeter Brune 
175e4094ef1SJacob Faibussowitsch   Developer Notes:
176f6dfbefdSBarry Smith   There should be a `SNESNCGSetType()`
177f6dfbefdSBarry Smith 
178*ceaaa498SBarry Smith .seealso: `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
1790a844d1aSPeter Brune @*/
180d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
181d71ae5a4SJacob Faibussowitsch {
1820a844d1aSPeter Brune   PetscFunctionBegin;
1830a844d1aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
184cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1860a844d1aSPeter Brune }
1870a844d1aSPeter Brune 
188d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
189d71ae5a4SJacob Faibussowitsch {
1900a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG *)snes->data;
1916e111a19SKarl Rupp 
1920a844d1aSPeter Brune   PetscFunctionBegin;
1930a844d1aSPeter Brune   ncg->type = btype;
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1950a844d1aSPeter Brune }
1960a844d1aSPeter Brune 
197fef7b6d8SPeter Brune /*
198fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
199fef7b6d8SPeter Brune 
2002fe279fdSBarry Smith   Input Parameter:
2012fe279fdSBarry Smith . snes - the `SNES` context
202fef7b6d8SPeter Brune 
203fef7b6d8SPeter Brune   Output Parameter:
204fef7b6d8SPeter Brune . outits - number of iterations until termination
205fef7b6d8SPeter Brune 
206fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
207fef7b6d8SPeter Brune */
208d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NCG(SNES snes)
209d71ae5a4SJacob Faibussowitsch {
210dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
21132cc618eSPeter Brune   Vec                  X, dX, lX, F, dXold;
212bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
21332cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
214fef7b6d8SPeter Brune   PetscInt             maxits, i;
215422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
216f1c6b773SPeter Brune   SNESLineSearch       linesearch;
217b7281c8aSPeter Brune   SNESConvergedReason  reason;
21888d374b2SPeter Brune 
219fef7b6d8SPeter Brune   PetscFunctionBegin;
2200b121fc5SBarry 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);
221c579b300SPatrick Farrell 
2229566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
223fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
224fef7b6d8SPeter Brune 
225fef7b6d8SPeter Brune   maxits = snes->max_its;        /* maximum number of iterations */
226fef7b6d8SPeter Brune   X      = snes->vec_sol;        /* X^n */
22732cc618eSPeter Brune   dXold  = snes->work[0];        /* The previous iterate of X */
22888d374b2SPeter Brune   dX     = snes->work[1];        /* the preconditioned direction */
229169e2be8SPeter Brune   lX     = snes->vec_sol_update; /* the conjugate direction */
230fef7b6d8SPeter Brune   F      = snes->vec_func;       /* residual vector */
231fef7b6d8SPeter Brune 
2329566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
2339e764e56SPeter Brune 
2349566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
235fef7b6d8SPeter Brune   snes->iter = 0;
236fef7b6d8SPeter Brune   snes->norm = 0.;
2379566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
238fef7b6d8SPeter Brune 
239bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
240a71552e2SPeter Brune 
241efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2429566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2439566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
244b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
245b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
2463ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
247b7281c8aSPeter Brune     }
2489566063dSJacob Faibussowitsch     PetscCall(VecCopy(dX, F));
2499566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
250a71552e2SPeter Brune   } else {
251e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2529566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2531aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
254c1c75074SPeter Brune 
255fef7b6d8SPeter Brune     /* convergence test */
2569566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
257422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2589566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, dX));
259a71552e2SPeter Brune   }
260efd4aadfSBarry Smith   if (snes->npc) {
261a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2629566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, dX));
2639566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
264b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
265b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2663ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
267b7281c8aSPeter Brune       }
268a71552e2SPeter Brune     }
269a71552e2SPeter Brune   }
2709566063dSJacob Faibussowitsch   PetscCall(VecCopy(dX, lX));
2719566063dSJacob Faibussowitsch   PetscCall(VecDot(dX, dX, &dXdotdX));
272a71552e2SPeter Brune 
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
274fef7b6d8SPeter Brune   snes->norm = fnorm;
2759566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2769566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
277fef7b6d8SPeter Brune 
278fef7b6d8SPeter Brune   /* test convergence */
2792d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
2802d157150SStefano Zampini   PetscCall(SNESMonitor(snes, 0, fnorm));
2813ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
282fef7b6d8SPeter Brune 
283fef7b6d8SPeter Brune   /* Call general purpose update function */
284dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
285fef7b6d8SPeter Brune 
286fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
287bbd5d0b3SPeter Brune 
28809c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
28929c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
29048a46eb9SPierre Jolivet     if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
2919566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
2929566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
2939566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
294422a814eSBarry Smith     if (lsresult) {
295fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
296fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2973ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
298fef7b6d8SPeter Brune       }
299fef7b6d8SPeter Brune     }
300e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
301fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3023ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
303fef7b6d8SPeter Brune     }
304fef7b6d8SPeter Brune     /* Monitor convergence */
3059566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
306169e2be8SPeter Brune     snes->iter  = i;
307fef7b6d8SPeter Brune     snes->norm  = fnorm;
308c1e67a49SFande Kong     snes->xnorm = xnorm;
309c1e67a49SFande Kong     snes->ynorm = ynorm;
3109566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3119566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
312302dbbaeSPeter Brune 
313fef7b6d8SPeter Brune     /* Test for convergence */
3142d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
3152d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
3163ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
317302dbbaeSPeter Brune 
318302dbbaeSPeter Brune     /* Call general purpose update function */
319dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
320efd4aadfSBarry Smith     if (snes->npc) {
321a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3229566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, dX));
3239566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
324b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
325b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
3263ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
327b7281c8aSPeter Brune         }
3289566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX, F));
329a71552e2SPeter Brune       } else {
3309566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, dX));
3319566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
332b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
333b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
3343ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
335b7281c8aSPeter Brune         }
336302dbbaeSPeter Brune       }
337a3685310SPeter Brune     } else {
3389566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, dX));
339302dbbaeSPeter Brune     }
340302dbbaeSPeter Brune 
34129c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3420a844d1aSPeter Brune     switch (ncg->type) {
3430a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
34432cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3459566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
34632cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / dXolddotdXold);
347dfb256c7SPeter Brune       break;
3480a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
34932cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3509566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
3519566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
3529566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3539566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
35432cc618eSPeter Brune       beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
355dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
356dfb256c7SPeter Brune       break;
3570a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
3589566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3599566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
3609566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3619566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3629566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3639566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
3649566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3659566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
36632cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
367dfb256c7SPeter Brune       break;
3680a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
3699566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3709566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3719566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3729566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3739566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3749566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
3752f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
376dfb256c7SPeter Brune       break;
3770a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
3789566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3799566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3809566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3819566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
3822f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
383dfb256c7SPeter Brune       break;
384dfb256c7SPeter Brune     }
38548a46eb9SPierre Jolivet     if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
3869566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
387fef7b6d8SPeter Brune   }
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
389fef7b6d8SPeter Brune }
390fef7b6d8SPeter Brune 
391fef7b6d8SPeter Brune /*MC
392f6dfbefdSBarry Smith   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems.
393fef7b6d8SPeter Brune 
394fef7b6d8SPeter Brune   Level: beginner
395fef7b6d8SPeter Brune 
396f6dfbefdSBarry Smith   Options Database Keys:
3974f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
39841a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
399b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
400fef7b6d8SPeter Brune 
40195452b02SPatrick Sanan    Notes:
40295452b02SPatrick Sanan    This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
403fef7b6d8SPeter Brune    gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
404fef7b6d8SPeter Brune    chooses the initial search direction as F(x) for the initial guess x.
405fef7b6d8SPeter Brune 
4062d547940SBarry Smith    Only supports left non-linear preconditioning.
4072d547940SBarry Smith 
408f6dfbefdSBarry Smith    Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then
409f6dfbefdSBarry Smith    `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR`
410b5badacbSBarry Smith 
4114f02bc6aSBarry Smith    References:
412606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4134f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
414cc3c3c0fSMatthew G. Knepley 
415f6dfbefdSBarry Smith .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
416fef7b6d8SPeter Brune M*/
417d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
418d71ae5a4SJacob Faibussowitsch {
419ea630c6eSPeter Brune   SNES_NCG *neP;
420fef7b6d8SPeter Brune 
421fef7b6d8SPeter Brune   PetscFunctionBegin;
422fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
423fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
424fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
425fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
426fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
427fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
428fef7b6d8SPeter Brune 
429fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
430efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
431efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
432fef7b6d8SPeter Brune 
4334fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4344fc747eaSLawrence Mitchell 
43588976e71SPeter Brune   if (!snes->tolerancesset) {
4360e444f03SPeter Brune     snes->max_funcs = 30000;
4370e444f03SPeter Brune     snes->max_its   = 10000;
438c522fa08SPeter Brune     snes->stol      = 1e-20;
43988976e71SPeter Brune   }
4400e444f03SPeter Brune 
4414dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
442fef7b6d8SPeter Brune   snes->data   = (void *)neP;
4430298fd71SBarry Smith   neP->monitor = NULL;
4440a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
4459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
4463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
447fef7b6d8SPeter Brune }
448