xref: /petsc/src/snes/impls/ncg/snesncg.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
49371c9d4SSatish Balay static PetscErrorCode SNESReset_NCG(SNES snes) {
5fef7b6d8SPeter Brune   PetscFunctionBegin;
6fef7b6d8SPeter Brune   PetscFunctionReturn(0);
7fef7b6d8SPeter Brune }
8fef7b6d8SPeter Brune 
9fef7b6d8SPeter Brune /*
10fef7b6d8SPeter Brune   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
11fef7b6d8SPeter Brune 
12fef7b6d8SPeter Brune   Input Parameter:
13fef7b6d8SPeter Brune . snes - the SNES context
14fef7b6d8SPeter Brune 
15fef7b6d8SPeter Brune   Application Interface Routine: SNESDestroy()
16fef7b6d8SPeter Brune */
179371c9d4SSatish Balay static PetscErrorCode SNESDestroy_NCG(SNES snes) {
18fef7b6d8SPeter Brune   PetscFunctionBegin;
192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
21fef7b6d8SPeter Brune   PetscFunctionReturn(0);
22fef7b6d8SPeter Brune }
23fef7b6d8SPeter Brune 
24fef7b6d8SPeter Brune /*
25fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
26fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
27fef7b6d8SPeter Brune 
28fef7b6d8SPeter Brune    Input Parameters:
29fef7b6d8SPeter Brune +  snes - the SNES context
30fef7b6d8SPeter Brune -  x - the solution vector
31fef7b6d8SPeter Brune 
32fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
33fef7b6d8SPeter Brune  */
34bbd5d0b3SPeter Brune 
359371c9d4SSatish Balay static PetscErrorCode SNESSetUp_NCG(SNES snes) {
36fef7b6d8SPeter Brune   PetscFunctionBegin;
379566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 2));
3808401ef6SPierre Jolivet   PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
396c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
40fef7b6d8SPeter Brune   PetscFunctionReturn(0);
41fef7b6d8SPeter Brune }
42fef7b6d8SPeter Brune 
439371c9d4SSatish Balay static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) {
44ea630c6eSPeter Brune   PetscScalar alpha, ptAp;
45bbd5d0b3SPeter Brune   Vec         X, Y, F, W;
46bbd5d0b3SPeter Brune   SNES        snes;
47bbd5d0b3SPeter Brune   PetscReal  *fnorm, *xnorm, *ynorm;
48ea630c6eSPeter Brune 
49ea630c6eSPeter Brune   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
51bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
52bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
53bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
54bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
55bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
56bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
57bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
58bbd5d0b3SPeter Brune 
5948a46eb9SPierre Jolivet   if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));
609105210eSPeter Brune 
61ea630c6eSPeter Brune   /*
62ea630c6eSPeter Brune 
63d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
64ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
65ea630c6eSPeter Brune    */
669566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
679566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &alpha));
689566063dSJacob Faibussowitsch   PetscCall(MatMult(snes->jacobian, Y, W));
699566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, W, &ptAp));
70ea630c6eSPeter Brune   alpha = alpha / ptAp;
719566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, -alpha, Y));
729566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, F));
73bbd5d0b3SPeter Brune 
749566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, fnorm));
759566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, xnorm));
769566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, ynorm));
77ea630c6eSPeter Brune   PetscFunctionReturn(0);
78ea630c6eSPeter Brune }
79ea630c6eSPeter Brune 
80b5badacbSBarry Smith /*MC
81*f6dfbefdSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG`
82b5badacbSBarry Smith 
83b5badacbSBarry 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.
84b5badacbSBarry 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.
85b5badacbSBarry Smith 
86*f6dfbefdSBarry Smith    Notes:
87*f6dfbefdSBarry Smith    This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
88b5badacbSBarry Smith 
89b5badacbSBarry Smith    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
90b5badacbSBarry Smith 
91b5badacbSBarry Smith    Level: advanced
92b5badacbSBarry Smith 
93db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
94b5badacbSBarry Smith M*/
95b5badacbSBarry Smith 
969371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) {
97bbd5d0b3SPeter Brune   PetscFunctionBegin;
98f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
990298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1000298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1010298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1020298fd71SBarry Smith   linesearch->ops->view           = NULL;
1030298fd71SBarry Smith   linesearch->ops->setup          = NULL;
104bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
105bbd5d0b3SPeter Brune }
106bbd5d0b3SPeter Brune 
10788d374b2SPeter Brune /*
108b5badacbSBarry Smith   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
109b5badacbSBarry Smith 
110b5badacbSBarry Smith   Input Parameter:
111b5badacbSBarry Smith . snes - the SNES context
112b5badacbSBarry Smith 
113b5badacbSBarry Smith   Application Interface Routine: SNESSetFromOptions()
114b5badacbSBarry Smith */
1159371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject) {
116b5badacbSBarry Smith   SNES_NCG      *ncg     = (SNES_NCG *)snes->data;
117b5badacbSBarry Smith   PetscBool      debug   = PETSC_FALSE;
118b5badacbSBarry Smith   SNESNCGType    ncgtype = ncg->type;
119b5badacbSBarry Smith   SNESLineSearch linesearch;
120b5badacbSBarry Smith 
121b5badacbSBarry Smith   PetscFunctionBegin;
122d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
1239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
124ad540459SPierre Jolivet   if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
1259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
1269566063dSJacob Faibussowitsch   PetscCall(SNESNCGSetType(snes, ncgtype));
127d0609cedSBarry Smith   PetscOptionsHeadEnd();
128b5badacbSBarry Smith   if (!snes->linesearch) {
1299566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
130ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
131b5badacbSBarry Smith       if (!snes->npc) {
1329566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
133b5badacbSBarry Smith       } else {
1349566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
135b5badacbSBarry Smith       }
136b5badacbSBarry Smith     }
137ec786807SJed Brown   }
138b5badacbSBarry Smith   PetscFunctionReturn(0);
139b5badacbSBarry Smith }
140b5badacbSBarry Smith 
141b5badacbSBarry Smith /*
142b5badacbSBarry Smith   SNESView_NCG - Prints info from the SNESNCG data structure.
143b5badacbSBarry Smith 
144b5badacbSBarry Smith   Input Parameters:
145b5badacbSBarry Smith + SNES - the SNES context
146b5badacbSBarry Smith - viewer - visualization context
147b5badacbSBarry Smith 
148b5badacbSBarry Smith   Application Interface Routine: SNESView()
149b5badacbSBarry Smith */
1509371c9d4SSatish Balay static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) {
151b5badacbSBarry Smith   SNES_NCG *ncg = (SNES_NCG *)snes->data;
152b5badacbSBarry Smith   PetscBool iascii;
153b5badacbSBarry Smith 
154b5badacbSBarry Smith   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
15648a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
157b5badacbSBarry Smith   PetscFunctionReturn(0);
158b5badacbSBarry Smith }
159b5badacbSBarry Smith 
160b5badacbSBarry Smith /*
16188d374b2SPeter Brune 
16288d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
16388d374b2SPeter Brune 
164b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
165b5badacbSBarry Smith 
166b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
167b5badacbSBarry Smith 
16888d374b2SPeter Brune  */
1699371c9d4SSatish Balay PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar *ytJtf) {
17088d374b2SPeter Brune   PetscScalar ftf, ftg, fty, h;
17167392de3SBarry Smith 
17288d374b2SPeter Brune   PetscFunctionBegin;
1739566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &ftf));
1749566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty));
17588d374b2SPeter Brune   h = 1e-5 * fty / fty;
1769566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
1779566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W, -h, Y)); /* this is arbitrary */
1789566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, W, G));
1799566063dSJacob Faibussowitsch   PetscCall(VecDot(G, F, &ftg));
18088d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
18188d374b2SPeter Brune   PetscFunctionReturn(0);
18288d374b2SPeter Brune }
18388d374b2SPeter Brune 
1840a844d1aSPeter Brune /*@
185*f6dfbefdSBarry Smith     SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`.
1860a844d1aSPeter Brune 
187*f6dfbefdSBarry Smith     Logically Collective on snes
1880a844d1aSPeter Brune 
1890a844d1aSPeter Brune     Input Parameters:
1900a844d1aSPeter Brune +   snes - the iterative context
1910a844d1aSPeter Brune -   btype - update type
1920a844d1aSPeter Brune 
1930a844d1aSPeter Brune     Options Database:
194b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
1950a844d1aSPeter Brune 
1960a844d1aSPeter Brune     Level: intermediate
1970a844d1aSPeter Brune 
198*f6dfbefdSBarry Smith     `SNESNCGType`s:
199*f6dfbefdSBarry Smith +   `SNES_NCG_FR` - Fletcher-Reeves update
200*f6dfbefdSBarry Smith .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update
201*f6dfbefdSBarry Smith .   `SNES_NCG_HS` - Hestenes-Steifel update
202*f6dfbefdSBarry Smith .   `SNES_NCG_DY` - Dai-Yuan update
203*f6dfbefdSBarry Smith -   `SNES_NCG_CD` - Conjugate Descent update
2040a844d1aSPeter Brune 
2050a844d1aSPeter Brune    Notes:
206*f6dfbefdSBarry Smith    `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.
207b5badacbSBarry Smith 
208b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
209*f6dfbefdSBarry Smith    that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`?
2100a844d1aSPeter Brune 
211*f6dfbefdSBarry Smith    Developer Note:
212*f6dfbefdSBarry Smith    There should be a `SNESNCGSetType()`
213*f6dfbefdSBarry Smith 
214*f6dfbefdSBarry Smith .seealso: `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`,  `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
2150a844d1aSPeter Brune @*/
2169371c9d4SSatish Balay PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) {
2170a844d1aSPeter Brune   PetscFunctionBegin;
2180a844d1aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
219cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
2200a844d1aSPeter Brune   PetscFunctionReturn(0);
2210a844d1aSPeter Brune }
2220a844d1aSPeter Brune 
2239371c9d4SSatish Balay static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) {
2240a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG *)snes->data;
2256e111a19SKarl Rupp 
2260a844d1aSPeter Brune   PetscFunctionBegin;
2270a844d1aSPeter Brune   ncg->type = btype;
2280a844d1aSPeter Brune   PetscFunctionReturn(0);
2290a844d1aSPeter Brune }
2300a844d1aSPeter Brune 
231fef7b6d8SPeter Brune /*
232fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
233fef7b6d8SPeter Brune 
234fef7b6d8SPeter Brune   Input Parameters:
235fef7b6d8SPeter Brune . snes - the SNES context
236fef7b6d8SPeter Brune 
237fef7b6d8SPeter Brune   Output Parameter:
238fef7b6d8SPeter Brune . outits - number of iterations until termination
239fef7b6d8SPeter Brune 
240fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
241fef7b6d8SPeter Brune */
2429371c9d4SSatish Balay static PetscErrorCode SNESSolve_NCG(SNES snes) {
243dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
24432cc618eSPeter Brune   Vec                  X, dX, lX, F, dXold;
245bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
24632cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
247fef7b6d8SPeter Brune   PetscInt             maxits, i;
248422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
249f1c6b773SPeter Brune   SNESLineSearch       linesearch;
250b7281c8aSPeter Brune   SNESConvergedReason  reason;
25188d374b2SPeter Brune 
252fef7b6d8SPeter Brune   PetscFunctionBegin;
2530b121fc5SBarry 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);
254c579b300SPatrick Farrell 
2559566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
256fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
257fef7b6d8SPeter Brune 
258fef7b6d8SPeter Brune   maxits = snes->max_its;        /* maximum number of iterations */
259fef7b6d8SPeter Brune   X      = snes->vec_sol;        /* X^n */
26032cc618eSPeter Brune   dXold  = snes->work[0];        /* The previous iterate of X */
26188d374b2SPeter Brune   dX     = snes->work[1];        /* the preconditioned direction */
262169e2be8SPeter Brune   lX     = snes->vec_sol_update; /* the conjugate direction */
263fef7b6d8SPeter Brune   F      = snes->vec_func;       /* residual vector */
264fef7b6d8SPeter Brune 
2659566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
2669e764e56SPeter Brune 
2679566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
268fef7b6d8SPeter Brune   snes->iter = 0;
269fef7b6d8SPeter Brune   snes->norm = 0.;
2709566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
271fef7b6d8SPeter Brune 
272bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
273a71552e2SPeter Brune 
274efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2759566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2769566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
277b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
278b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
279b7281c8aSPeter Brune       PetscFunctionReturn(0);
280b7281c8aSPeter Brune     }
2819566063dSJacob Faibussowitsch     PetscCall(VecCopy(dX, F));
2829566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
283a71552e2SPeter Brune   } else {
284e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2859566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2861aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
287c1c75074SPeter Brune 
288fef7b6d8SPeter Brune     /* convergence test */
2899566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
290422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2919566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, dX));
292a71552e2SPeter Brune   }
293efd4aadfSBarry Smith   if (snes->npc) {
294a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2959566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, dX));
2969566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
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     }
302a71552e2SPeter Brune   }
3039566063dSJacob Faibussowitsch   PetscCall(VecCopy(dX, lX));
3049566063dSJacob Faibussowitsch   PetscCall(VecDot(dX, dX, &dXdotdX));
305a71552e2SPeter Brune 
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
307fef7b6d8SPeter Brune   snes->norm = fnorm;
3089566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3099566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3109566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
311fef7b6d8SPeter Brune 
312fef7b6d8SPeter Brune   /* test convergence */
313dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
314fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
315fef7b6d8SPeter Brune 
316fef7b6d8SPeter Brune   /* Call general purpose update function */
317dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
318fef7b6d8SPeter Brune 
319fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
320bbd5d0b3SPeter Brune 
32109c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
32229c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
32348a46eb9SPierre Jolivet     if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
3249566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
3259566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
3269566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
327422a814eSBarry Smith     if (lsresult) {
328fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
329fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3305d10c001SPeter Brune         PetscFunctionReturn(0);
331fef7b6d8SPeter Brune       }
332fef7b6d8SPeter Brune     }
333e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
334fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3355d10c001SPeter Brune       PetscFunctionReturn(0);
336fef7b6d8SPeter Brune     }
337fef7b6d8SPeter Brune     /* Monitor convergence */
3389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
339169e2be8SPeter Brune     snes->iter  = i;
340fef7b6d8SPeter Brune     snes->norm  = fnorm;
341c1e67a49SFande Kong     snes->xnorm = xnorm;
342c1e67a49SFande Kong     snes->ynorm = ynorm;
3439566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3449566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
3459566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
346302dbbaeSPeter Brune 
347fef7b6d8SPeter Brune     /* Test for convergence */
348dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
3495d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
350302dbbaeSPeter Brune 
351302dbbaeSPeter Brune     /* Call general purpose update function */
352dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
353efd4aadfSBarry Smith     if (snes->npc) {
354a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3559566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, dX));
3569566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
357b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
358b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
359b7281c8aSPeter Brune           PetscFunctionReturn(0);
360b7281c8aSPeter Brune         }
3619566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX, F));
362a71552e2SPeter Brune       } else {
3639566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, dX));
3649566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
365b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
366b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
367b7281c8aSPeter Brune           PetscFunctionReturn(0);
368b7281c8aSPeter Brune         }
369302dbbaeSPeter Brune       }
370a3685310SPeter Brune     } else {
3719566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, dX));
372302dbbaeSPeter Brune     }
373302dbbaeSPeter Brune 
37429c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3750a844d1aSPeter Brune     switch (ncg->type) {
3760a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
37732cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3789566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
37932cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / dXolddotdXold);
380dfb256c7SPeter Brune       break;
3810a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
38232cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3839566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
3849566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
3859566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3869566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
38732cc618eSPeter Brune       beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
388dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
389dfb256c7SPeter Brune       break;
3900a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
3919566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3929566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
3939566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3949566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3959566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3969566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
3979566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3989566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
39932cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
400dfb256c7SPeter Brune       break;
4010a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
4029566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4039566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4049566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4059566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4069566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4079566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4082f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
409dfb256c7SPeter Brune       break;
4100a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
4119566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4129566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4139566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4149566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4152f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
416dfb256c7SPeter Brune       break;
417dfb256c7SPeter Brune     }
41848a46eb9SPierre Jolivet     if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
4199566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
420fef7b6d8SPeter Brune   }
42163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
422fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
423fef7b6d8SPeter Brune   PetscFunctionReturn(0);
424fef7b6d8SPeter Brune }
425fef7b6d8SPeter Brune 
426fef7b6d8SPeter Brune /*MC
427*f6dfbefdSBarry Smith   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems.
428fef7b6d8SPeter Brune 
429fef7b6d8SPeter Brune   Level: beginner
430fef7b6d8SPeter Brune 
431*f6dfbefdSBarry Smith   Options Database Keys:
4324f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
43341a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
434b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
435fef7b6d8SPeter Brune 
43695452b02SPatrick Sanan    Notes:
43795452b02SPatrick Sanan    This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
438fef7b6d8SPeter Brune    gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
439fef7b6d8SPeter Brune    chooses the initial search direction as F(x) for the initial guess x.
440fef7b6d8SPeter Brune 
4412d547940SBarry Smith    Only supports left non-linear preconditioning.
4422d547940SBarry Smith 
443*f6dfbefdSBarry Smith    Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then
444*f6dfbefdSBarry Smith    `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR`
445b5badacbSBarry Smith 
4464f02bc6aSBarry Smith    References:
447606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4484f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
449cc3c3c0fSMatthew G. Knepley 
450*f6dfbefdSBarry Smith .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
451fef7b6d8SPeter Brune M*/
4529371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) {
453ea630c6eSPeter Brune   SNES_NCG *neP;
454fef7b6d8SPeter Brune 
455fef7b6d8SPeter Brune   PetscFunctionBegin;
456fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
457fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
458fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
459fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
460fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
461fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
462fef7b6d8SPeter Brune 
463fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
464efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
465efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
466fef7b6d8SPeter Brune 
4674fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4684fc747eaSLawrence Mitchell 
46988976e71SPeter Brune   if (!snes->tolerancesset) {
4700e444f03SPeter Brune     snes->max_funcs = 30000;
4710e444f03SPeter Brune     snes->max_its   = 10000;
472c522fa08SPeter Brune     snes->stol      = 1e-20;
47388976e71SPeter Brune   }
4740e444f03SPeter Brune 
4759566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &neP));
476fef7b6d8SPeter Brune   snes->data   = (void *)neP;
4770298fd71SBarry Smith   neP->monitor = NULL;
4780a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
4799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
480fef7b6d8SPeter Brune   PetscFunctionReturn(0);
481fef7b6d8SPeter Brune }
482