xref: /petsc/src/snes/impls/ncg/snesncg.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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;
7fef7b6d8SPeter Brune   PetscFunctionReturn(0);
8fef7b6d8SPeter Brune }
9fef7b6d8SPeter Brune 
10fef7b6d8SPeter Brune /*
11fef7b6d8SPeter Brune   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
12fef7b6d8SPeter Brune 
13fef7b6d8SPeter Brune   Input Parameter:
14fef7b6d8SPeter Brune . snes - the SNES context
15fef7b6d8SPeter Brune 
16fef7b6d8SPeter Brune   Application Interface Routine: SNESDestroy()
17fef7b6d8SPeter Brune */
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NCG(SNES snes)
19d71ae5a4SJacob Faibussowitsch {
20fef7b6d8SPeter Brune   PetscFunctionBegin;
212e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
23fef7b6d8SPeter Brune   PetscFunctionReturn(0);
24fef7b6d8SPeter Brune }
25fef7b6d8SPeter Brune 
26fef7b6d8SPeter Brune /*
27fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
28fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
29fef7b6d8SPeter Brune 
30fef7b6d8SPeter Brune    Input Parameters:
31fef7b6d8SPeter Brune +  snes - the SNES context
32fef7b6d8SPeter Brune -  x - the solution vector
33fef7b6d8SPeter Brune 
34fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
35fef7b6d8SPeter Brune  */
36bbd5d0b3SPeter Brune 
37d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NCG(SNES snes)
38d71ae5a4SJacob Faibussowitsch {
39fef7b6d8SPeter Brune   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 2));
4108401ef6SPierre Jolivet   PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
426c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
43fef7b6d8SPeter Brune   PetscFunctionReturn(0);
44fef7b6d8SPeter Brune }
45fef7b6d8SPeter Brune 
46d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
47d71ae5a4SJacob Faibussowitsch {
48ea630c6eSPeter Brune   PetscScalar alpha, ptAp;
49bbd5d0b3SPeter Brune   Vec         X, Y, F, W;
50bbd5d0b3SPeter Brune   SNES        snes;
51bbd5d0b3SPeter Brune   PetscReal  *fnorm, *xnorm, *ynorm;
52ea630c6eSPeter Brune 
53ea630c6eSPeter Brune   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
55bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
56bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
57bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
58bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
59bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
60bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
61bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
62bbd5d0b3SPeter Brune 
6348a46eb9SPierre Jolivet   if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));
649105210eSPeter Brune 
65ea630c6eSPeter Brune   /*
66ea630c6eSPeter Brune 
67d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
68ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
69ea630c6eSPeter Brune    */
709566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
719566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &alpha));
729566063dSJacob Faibussowitsch   PetscCall(MatMult(snes->jacobian, Y, W));
739566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, W, &ptAp));
74ea630c6eSPeter Brune   alpha = alpha / ptAp;
759566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, -alpha, Y));
769566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, F));
77bbd5d0b3SPeter Brune 
789566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, fnorm));
799566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, xnorm));
809566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y, NORM_2, ynorm));
81ea630c6eSPeter Brune   PetscFunctionReturn(0);
82ea630c6eSPeter Brune }
83ea630c6eSPeter Brune 
84b5badacbSBarry Smith /*MC
85f6dfbefdSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG`
86b5badacbSBarry Smith 
87b5badacbSBarry 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.
88b5badacbSBarry 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.
89b5badacbSBarry Smith 
90f6dfbefdSBarry Smith    Notes:
91f6dfbefdSBarry Smith    This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
92b5badacbSBarry Smith 
93b5badacbSBarry Smith    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
94b5badacbSBarry Smith 
95b5badacbSBarry Smith    Level: advanced
96b5badacbSBarry Smith 
97db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
98b5badacbSBarry Smith M*/
99b5badacbSBarry Smith 
100d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
101d71ae5a4SJacob Faibussowitsch {
102bbd5d0b3SPeter Brune   PetscFunctionBegin;
103f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
1040298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1050298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1060298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1070298fd71SBarry Smith   linesearch->ops->view           = NULL;
1080298fd71SBarry Smith   linesearch->ops->setup          = NULL;
109bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
110bbd5d0b3SPeter Brune }
111bbd5d0b3SPeter Brune 
11288d374b2SPeter Brune /*
113b5badacbSBarry Smith   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
114b5badacbSBarry Smith 
115b5badacbSBarry Smith   Input Parameter:
116b5badacbSBarry Smith . snes - the SNES context
117b5badacbSBarry Smith 
118b5badacbSBarry Smith   Application Interface Routine: SNESSetFromOptions()
119b5badacbSBarry Smith */
120d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject)
121d71ae5a4SJacob Faibussowitsch {
122b5badacbSBarry Smith   SNES_NCG      *ncg     = (SNES_NCG *)snes->data;
123b5badacbSBarry Smith   PetscBool      debug   = PETSC_FALSE;
124b5badacbSBarry Smith   SNESNCGType    ncgtype = ncg->type;
125b5badacbSBarry Smith   SNESLineSearch linesearch;
126b5badacbSBarry Smith 
127b5badacbSBarry Smith   PetscFunctionBegin;
128d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
1299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
130ad540459SPierre Jolivet   if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
1319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
1329566063dSJacob Faibussowitsch   PetscCall(SNESNCGSetType(snes, ncgtype));
133d0609cedSBarry Smith   PetscOptionsHeadEnd();
134b5badacbSBarry Smith   if (!snes->linesearch) {
1359566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
136ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
137b5badacbSBarry Smith       if (!snes->npc) {
1389566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
139b5badacbSBarry Smith       } else {
1409566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
141b5badacbSBarry Smith       }
142b5badacbSBarry Smith     }
143ec786807SJed Brown   }
144b5badacbSBarry Smith   PetscFunctionReturn(0);
145b5badacbSBarry Smith }
146b5badacbSBarry Smith 
147b5badacbSBarry Smith /*
148b5badacbSBarry Smith   SNESView_NCG - Prints info from the SNESNCG data structure.
149b5badacbSBarry Smith 
150b5badacbSBarry Smith   Input Parameters:
151b5badacbSBarry Smith + SNES - the SNES context
152b5badacbSBarry Smith - viewer - visualization context
153b5badacbSBarry Smith 
154b5badacbSBarry Smith   Application Interface Routine: SNESView()
155b5badacbSBarry Smith */
156d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
157d71ae5a4SJacob Faibussowitsch {
158b5badacbSBarry Smith   SNES_NCG *ncg = (SNES_NCG *)snes->data;
159b5badacbSBarry Smith   PetscBool iascii;
160b5badacbSBarry Smith 
161b5badacbSBarry Smith   PetscFunctionBegin;
1629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
16348a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
164b5badacbSBarry Smith   PetscFunctionReturn(0);
165b5badacbSBarry Smith }
166b5badacbSBarry Smith 
167b5badacbSBarry Smith /*
16888d374b2SPeter Brune 
16988d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
17088d374b2SPeter Brune 
171b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
172b5badacbSBarry Smith 
173b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
174b5badacbSBarry Smith 
17588d374b2SPeter Brune  */
176d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar *ytJtf)
177d71ae5a4SJacob Faibussowitsch {
17888d374b2SPeter Brune   PetscScalar ftf, ftg, fty, h;
17967392de3SBarry Smith 
18088d374b2SPeter Brune   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &ftf));
1829566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty));
18388d374b2SPeter Brune   h = 1e-5 * fty / fty;
1849566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
1859566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W, -h, Y)); /* this is arbitrary */
1869566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, W, G));
1879566063dSJacob Faibussowitsch   PetscCall(VecDot(G, F, &ftg));
18888d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
18988d374b2SPeter Brune   PetscFunctionReturn(0);
19088d374b2SPeter Brune }
19188d374b2SPeter Brune 
1920a844d1aSPeter Brune /*@
193f6dfbefdSBarry Smith     SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`.
1940a844d1aSPeter Brune 
195*c3339decSBarry Smith     Logically Collective
1960a844d1aSPeter Brune 
1970a844d1aSPeter Brune     Input Parameters:
1980a844d1aSPeter Brune +   snes - the iterative context
1990a844d1aSPeter Brune -   btype - update type
2000a844d1aSPeter Brune 
2013c7db156SBarry Smith     Options Database Key:
202b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
2030a844d1aSPeter Brune 
2040a844d1aSPeter Brune     Level: intermediate
2050a844d1aSPeter Brune 
206f6dfbefdSBarry Smith     `SNESNCGType`s:
207f6dfbefdSBarry Smith +   `SNES_NCG_FR` - Fletcher-Reeves update
208f6dfbefdSBarry Smith .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update
209f6dfbefdSBarry Smith .   `SNES_NCG_HS` - Hestenes-Steifel update
210f6dfbefdSBarry Smith .   `SNES_NCG_DY` - Dai-Yuan update
211f6dfbefdSBarry Smith -   `SNES_NCG_CD` - Conjugate Descent update
2120a844d1aSPeter Brune 
2130a844d1aSPeter Brune    Notes:
214f6dfbefdSBarry Smith    `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.
215b5badacbSBarry Smith 
216b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
217f6dfbefdSBarry Smith    that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`?
2180a844d1aSPeter Brune 
219f6dfbefdSBarry Smith    Developer Note:
220f6dfbefdSBarry Smith    There should be a `SNESNCGSetType()`
221f6dfbefdSBarry Smith 
222f6dfbefdSBarry Smith .seealso: `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
2230a844d1aSPeter Brune @*/
224d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
225d71ae5a4SJacob Faibussowitsch {
2260a844d1aSPeter Brune   PetscFunctionBegin;
2270a844d1aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
228cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
2290a844d1aSPeter Brune   PetscFunctionReturn(0);
2300a844d1aSPeter Brune }
2310a844d1aSPeter Brune 
232d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
233d71ae5a4SJacob Faibussowitsch {
2340a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG *)snes->data;
2356e111a19SKarl Rupp 
2360a844d1aSPeter Brune   PetscFunctionBegin;
2370a844d1aSPeter Brune   ncg->type = btype;
2380a844d1aSPeter Brune   PetscFunctionReturn(0);
2390a844d1aSPeter Brune }
2400a844d1aSPeter Brune 
241fef7b6d8SPeter Brune /*
242fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
243fef7b6d8SPeter Brune 
244fef7b6d8SPeter Brune   Input Parameters:
245fef7b6d8SPeter Brune . snes - the SNES context
246fef7b6d8SPeter Brune 
247fef7b6d8SPeter Brune   Output Parameter:
248fef7b6d8SPeter Brune . outits - number of iterations until termination
249fef7b6d8SPeter Brune 
250fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
251fef7b6d8SPeter Brune */
252d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NCG(SNES snes)
253d71ae5a4SJacob Faibussowitsch {
254dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
25532cc618eSPeter Brune   Vec                  X, dX, lX, F, dXold;
256bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
25732cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
258fef7b6d8SPeter Brune   PetscInt             maxits, i;
259422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
260f1c6b773SPeter Brune   SNESLineSearch       linesearch;
261b7281c8aSPeter Brune   SNESConvergedReason  reason;
26288d374b2SPeter Brune 
263fef7b6d8SPeter Brune   PetscFunctionBegin;
2640b121fc5SBarry 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);
265c579b300SPatrick Farrell 
2669566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
267fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
268fef7b6d8SPeter Brune 
269fef7b6d8SPeter Brune   maxits = snes->max_its;        /* maximum number of iterations */
270fef7b6d8SPeter Brune   X      = snes->vec_sol;        /* X^n */
27132cc618eSPeter Brune   dXold  = snes->work[0];        /* The previous iterate of X */
27288d374b2SPeter Brune   dX     = snes->work[1];        /* the preconditioned direction */
273169e2be8SPeter Brune   lX     = snes->vec_sol_update; /* the conjugate direction */
274fef7b6d8SPeter Brune   F      = snes->vec_func;       /* residual vector */
275fef7b6d8SPeter Brune 
2769566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
2779e764e56SPeter Brune 
2789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
279fef7b6d8SPeter Brune   snes->iter = 0;
280fef7b6d8SPeter Brune   snes->norm = 0.;
2819566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
282fef7b6d8SPeter Brune 
283bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
284a71552e2SPeter Brune 
285efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2869566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2879566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
288b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
289b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
290b7281c8aSPeter Brune       PetscFunctionReturn(0);
291b7281c8aSPeter Brune     }
2929566063dSJacob Faibussowitsch     PetscCall(VecCopy(dX, F));
2939566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
294a71552e2SPeter Brune   } else {
295e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2969566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2971aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
298c1c75074SPeter Brune 
299fef7b6d8SPeter Brune     /* convergence test */
3009566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
301422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
3029566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, dX));
303a71552e2SPeter Brune   }
304efd4aadfSBarry Smith   if (snes->npc) {
305a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
3069566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, dX));
3079566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
308b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
309b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
310b7281c8aSPeter Brune         PetscFunctionReturn(0);
311b7281c8aSPeter Brune       }
312a71552e2SPeter Brune     }
313a71552e2SPeter Brune   }
3149566063dSJacob Faibussowitsch   PetscCall(VecCopy(dX, lX));
3159566063dSJacob Faibussowitsch   PetscCall(VecDot(dX, dX, &dXdotdX));
316a71552e2SPeter Brune 
3179566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
318fef7b6d8SPeter Brune   snes->norm = fnorm;
3199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3209566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3219566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
322fef7b6d8SPeter Brune 
323fef7b6d8SPeter Brune   /* test convergence */
324dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
325fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
326fef7b6d8SPeter Brune 
327fef7b6d8SPeter Brune   /* Call general purpose update function */
328dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
329fef7b6d8SPeter Brune 
330fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
331bbd5d0b3SPeter Brune 
33209c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
33329c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
33448a46eb9SPierre Jolivet     if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
3359566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
3369566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
3379566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
338422a814eSBarry Smith     if (lsresult) {
339fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
340fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3415d10c001SPeter Brune         PetscFunctionReturn(0);
342fef7b6d8SPeter Brune       }
343fef7b6d8SPeter Brune     }
344e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
345fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3465d10c001SPeter Brune       PetscFunctionReturn(0);
347fef7b6d8SPeter Brune     }
348fef7b6d8SPeter Brune     /* Monitor convergence */
3499566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
350169e2be8SPeter Brune     snes->iter  = i;
351fef7b6d8SPeter Brune     snes->norm  = fnorm;
352c1e67a49SFande Kong     snes->xnorm = xnorm;
353c1e67a49SFande Kong     snes->ynorm = ynorm;
3549566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3559566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
3569566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
357302dbbaeSPeter Brune 
358fef7b6d8SPeter Brune     /* Test for convergence */
359dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
3605d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
361302dbbaeSPeter Brune 
362302dbbaeSPeter Brune     /* Call general purpose update function */
363dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
364efd4aadfSBarry Smith     if (snes->npc) {
365a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3669566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, dX));
3679566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
368b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
369b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
370b7281c8aSPeter Brune           PetscFunctionReturn(0);
371b7281c8aSPeter Brune         }
3729566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX, F));
373a71552e2SPeter Brune       } else {
3749566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, dX));
3759566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
376b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
377b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
378b7281c8aSPeter Brune           PetscFunctionReturn(0);
379b7281c8aSPeter Brune         }
380302dbbaeSPeter Brune       }
381a3685310SPeter Brune     } else {
3829566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, dX));
383302dbbaeSPeter Brune     }
384302dbbaeSPeter Brune 
38529c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3860a844d1aSPeter Brune     switch (ncg->type) {
3870a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
38832cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3899566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
39032cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / dXolddotdXold);
391dfb256c7SPeter Brune       break;
3920a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
39332cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3949566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
3959566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
3969566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3979566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
39832cc618eSPeter Brune       beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
399dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
400dfb256c7SPeter Brune       break;
4010a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
4029566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4039566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
4049566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4059566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4069566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4079566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
4089566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4099566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
41032cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
411dfb256c7SPeter Brune       break;
4120a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
4139566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4149566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4159566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4169566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4179566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4189566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4192f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
420dfb256c7SPeter Brune       break;
4210a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
4229566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4239566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4249566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4259566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4262f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
427dfb256c7SPeter Brune       break;
428dfb256c7SPeter Brune     }
42948a46eb9SPierre Jolivet     if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
4309566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
431fef7b6d8SPeter Brune   }
43263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
433fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
434fef7b6d8SPeter Brune   PetscFunctionReturn(0);
435fef7b6d8SPeter Brune }
436fef7b6d8SPeter Brune 
437fef7b6d8SPeter Brune /*MC
438f6dfbefdSBarry Smith   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems.
439fef7b6d8SPeter Brune 
440fef7b6d8SPeter Brune   Level: beginner
441fef7b6d8SPeter Brune 
442f6dfbefdSBarry Smith   Options Database Keys:
4434f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
44441a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
445b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
446fef7b6d8SPeter Brune 
44795452b02SPatrick Sanan    Notes:
44895452b02SPatrick Sanan    This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
449fef7b6d8SPeter Brune    gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
450fef7b6d8SPeter Brune    chooses the initial search direction as F(x) for the initial guess x.
451fef7b6d8SPeter Brune 
4522d547940SBarry Smith    Only supports left non-linear preconditioning.
4532d547940SBarry Smith 
454f6dfbefdSBarry Smith    Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then
455f6dfbefdSBarry Smith    `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR`
456b5badacbSBarry Smith 
4574f02bc6aSBarry Smith    References:
458606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4594f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
460cc3c3c0fSMatthew G. Knepley 
461f6dfbefdSBarry Smith .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
462fef7b6d8SPeter Brune M*/
463d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
464d71ae5a4SJacob Faibussowitsch {
465ea630c6eSPeter Brune   SNES_NCG *neP;
466fef7b6d8SPeter Brune 
467fef7b6d8SPeter Brune   PetscFunctionBegin;
468fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
469fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
470fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
471fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
472fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
473fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
474fef7b6d8SPeter Brune 
475fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
476efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
477efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
478fef7b6d8SPeter Brune 
4794fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4804fc747eaSLawrence Mitchell 
48188976e71SPeter Brune   if (!snes->tolerancesset) {
4820e444f03SPeter Brune     snes->max_funcs = 30000;
4830e444f03SPeter Brune     snes->max_its   = 10000;
484c522fa08SPeter Brune     snes->stol      = 1e-20;
48588976e71SPeter Brune   }
4860e444f03SPeter Brune 
4874dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
488fef7b6d8SPeter Brune   snes->data   = (void *)neP;
4890298fd71SBarry Smith   neP->monitor = NULL;
4900a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
4919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
492fef7b6d8SPeter Brune   PetscFunctionReturn(0);
493fef7b6d8SPeter Brune }
494