xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
59*48a46eb9SPierre 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
81b5badacbSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for 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 
86b5badacbSBarry Smith    Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
87b5badacbSBarry Smith 
88b5badacbSBarry Smith    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
89b5badacbSBarry Smith 
90b5badacbSBarry Smith    Level: advanced
91b5badacbSBarry Smith 
92db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
93b5badacbSBarry Smith M*/
94b5badacbSBarry Smith 
959371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) {
96bbd5d0b3SPeter Brune   PetscFunctionBegin;
97f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
980298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
990298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1000298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1010298fd71SBarry Smith   linesearch->ops->view           = NULL;
1020298fd71SBarry Smith   linesearch->ops->setup          = NULL;
103bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
104bbd5d0b3SPeter Brune }
105bbd5d0b3SPeter Brune 
10688d374b2SPeter Brune /*
107b5badacbSBarry Smith   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
108b5badacbSBarry Smith 
109b5badacbSBarry Smith   Input Parameter:
110b5badacbSBarry Smith . snes - the SNES context
111b5badacbSBarry Smith 
112b5badacbSBarry Smith   Application Interface Routine: SNESSetFromOptions()
113b5badacbSBarry Smith */
1149371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject) {
115b5badacbSBarry Smith   SNES_NCG      *ncg     = (SNES_NCG *)snes->data;
116b5badacbSBarry Smith   PetscBool      debug   = PETSC_FALSE;
117b5badacbSBarry Smith   SNESNCGType    ncgtype = ncg->type;
118b5badacbSBarry Smith   SNESLineSearch linesearch;
119b5badacbSBarry Smith 
120b5badacbSBarry Smith   PetscFunctionBegin;
121d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
1229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
1239371c9d4SSatish Balay   if (debug) { ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); }
1249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
1259566063dSJacob Faibussowitsch   PetscCall(SNESNCGSetType(snes, ncgtype));
126d0609cedSBarry Smith   PetscOptionsHeadEnd();
127b5badacbSBarry Smith   if (!snes->linesearch) {
1289566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
129ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
130b5badacbSBarry Smith       if (!snes->npc) {
1319566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
132b5badacbSBarry Smith       } else {
1339566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
134b5badacbSBarry Smith       }
135b5badacbSBarry Smith     }
136ec786807SJed Brown   }
137b5badacbSBarry Smith   PetscFunctionReturn(0);
138b5badacbSBarry Smith }
139b5badacbSBarry Smith 
140b5badacbSBarry Smith /*
141b5badacbSBarry Smith   SNESView_NCG - Prints info from the SNESNCG data structure.
142b5badacbSBarry Smith 
143b5badacbSBarry Smith   Input Parameters:
144b5badacbSBarry Smith + SNES - the SNES context
145b5badacbSBarry Smith - viewer - visualization context
146b5badacbSBarry Smith 
147b5badacbSBarry Smith   Application Interface Routine: SNESView()
148b5badacbSBarry Smith */
1499371c9d4SSatish Balay static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) {
150b5badacbSBarry Smith   SNES_NCG *ncg = (SNES_NCG *)snes->data;
151b5badacbSBarry Smith   PetscBool iascii;
152b5badacbSBarry Smith 
153b5badacbSBarry Smith   PetscFunctionBegin;
1549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
155*48a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
156b5badacbSBarry Smith   PetscFunctionReturn(0);
157b5badacbSBarry Smith }
158b5badacbSBarry Smith 
159b5badacbSBarry Smith /*
16088d374b2SPeter Brune 
16188d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
16288d374b2SPeter Brune 
163b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
164b5badacbSBarry Smith 
165b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
166b5badacbSBarry Smith 
16788d374b2SPeter Brune  */
1689371c9d4SSatish Balay PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar *ytJtf) {
16988d374b2SPeter Brune   PetscScalar ftf, ftg, fty, h;
17067392de3SBarry Smith 
17188d374b2SPeter Brune   PetscFunctionBegin;
1729566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &ftf));
1739566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty));
17488d374b2SPeter Brune   h = 1e-5 * fty / fty;
1759566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
1769566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W, -h, Y)); /* this is arbitrary */
1779566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, W, G));
1789566063dSJacob Faibussowitsch   PetscCall(VecDot(G, F, &ftg));
17988d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
18088d374b2SPeter Brune   PetscFunctionReturn(0);
18188d374b2SPeter Brune }
18288d374b2SPeter Brune 
1830a844d1aSPeter Brune /*@
1840a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
1850a844d1aSPeter Brune 
1860a844d1aSPeter Brune     Logically Collective on SNES
1870a844d1aSPeter Brune 
1880a844d1aSPeter Brune     Input Parameters:
1890a844d1aSPeter Brune +   snes - the iterative context
1900a844d1aSPeter Brune -   btype - update type
1910a844d1aSPeter Brune 
1920a844d1aSPeter Brune     Options Database:
193b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
1940a844d1aSPeter Brune 
1950a844d1aSPeter Brune     Level: intermediate
1960a844d1aSPeter Brune 
1970a844d1aSPeter Brune     SNESNCGSelectTypes:
1980a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
1990a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2000a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2010a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2020a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2030a844d1aSPeter Brune 
2040a844d1aSPeter Brune    Notes:
205b5badacbSBarry Smith    SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.
206b5badacbSBarry Smith 
207b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
208b5badacbSBarry Smith    that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?
2090a844d1aSPeter Brune 
2100a844d1aSPeter Brune @*/
2119371c9d4SSatish Balay PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) {
2120a844d1aSPeter Brune   PetscFunctionBegin;
2130a844d1aSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
214cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
2150a844d1aSPeter Brune   PetscFunctionReturn(0);
2160a844d1aSPeter Brune }
2170a844d1aSPeter Brune 
2189371c9d4SSatish Balay static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) {
2190a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG *)snes->data;
2206e111a19SKarl Rupp 
2210a844d1aSPeter Brune   PetscFunctionBegin;
2220a844d1aSPeter Brune   ncg->type = btype;
2230a844d1aSPeter Brune   PetscFunctionReturn(0);
2240a844d1aSPeter Brune }
2250a844d1aSPeter Brune 
226fef7b6d8SPeter Brune /*
227fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
228fef7b6d8SPeter Brune 
229fef7b6d8SPeter Brune   Input Parameters:
230fef7b6d8SPeter Brune . snes - the SNES context
231fef7b6d8SPeter Brune 
232fef7b6d8SPeter Brune   Output Parameter:
233fef7b6d8SPeter Brune . outits - number of iterations until termination
234fef7b6d8SPeter Brune 
235fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
236fef7b6d8SPeter Brune */
2379371c9d4SSatish Balay static PetscErrorCode SNESSolve_NCG(SNES snes) {
238dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
23932cc618eSPeter Brune   Vec                  X, dX, lX, F, dXold;
240bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
24132cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
242fef7b6d8SPeter Brune   PetscInt             maxits, i;
243422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
244f1c6b773SPeter Brune   SNESLineSearch       linesearch;
245b7281c8aSPeter Brune   SNESConvergedReason  reason;
24688d374b2SPeter Brune 
247fef7b6d8SPeter Brune   PetscFunctionBegin;
2480b121fc5SBarry 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);
249c579b300SPatrick Farrell 
2509566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
251fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
252fef7b6d8SPeter Brune 
253fef7b6d8SPeter Brune   maxits = snes->max_its;        /* maximum number of iterations */
254fef7b6d8SPeter Brune   X      = snes->vec_sol;        /* X^n */
25532cc618eSPeter Brune   dXold  = snes->work[0];        /* The previous iterate of X */
25688d374b2SPeter Brune   dX     = snes->work[1];        /* the preconditioned direction */
257169e2be8SPeter Brune   lX     = snes->vec_sol_update; /* the conjugate direction */
258fef7b6d8SPeter Brune   F      = snes->vec_func;       /* residual vector */
259fef7b6d8SPeter Brune 
2609566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
2619e764e56SPeter Brune 
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
263fef7b6d8SPeter Brune   snes->iter = 0;
264fef7b6d8SPeter Brune   snes->norm = 0.;
2659566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
266fef7b6d8SPeter Brune 
267bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
268a71552e2SPeter Brune 
269efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2709566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2719566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
272b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
273b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
274b7281c8aSPeter Brune       PetscFunctionReturn(0);
275b7281c8aSPeter Brune     }
2769566063dSJacob Faibussowitsch     PetscCall(VecCopy(dX, F));
2779566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
278a71552e2SPeter Brune   } else {
279e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2809566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2811aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
282c1c75074SPeter Brune 
283fef7b6d8SPeter Brune     /* convergence test */
2849566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
285422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2869566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, dX));
287a71552e2SPeter Brune   }
288efd4aadfSBarry Smith   if (snes->npc) {
289a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2909566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, X, F, dX));
2919566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
292b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
293b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
294b7281c8aSPeter Brune         PetscFunctionReturn(0);
295b7281c8aSPeter Brune       }
296a71552e2SPeter Brune     }
297a71552e2SPeter Brune   }
2989566063dSJacob Faibussowitsch   PetscCall(VecCopy(dX, lX));
2999566063dSJacob Faibussowitsch   PetscCall(VecDot(dX, dX, &dXdotdX));
300a71552e2SPeter Brune 
3019566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
302fef7b6d8SPeter Brune   snes->norm = fnorm;
3039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3049566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3059566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
306fef7b6d8SPeter Brune 
307fef7b6d8SPeter Brune   /* test convergence */
308dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
309fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
310fef7b6d8SPeter Brune 
311fef7b6d8SPeter Brune   /* Call general purpose update function */
312dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
313fef7b6d8SPeter Brune 
314fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
315bbd5d0b3SPeter Brune 
31609c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
31729c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
318*48a46eb9SPierre Jolivet     if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
3199566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
3209566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
3219566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
322422a814eSBarry Smith     if (lsresult) {
323fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
324fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3255d10c001SPeter Brune         PetscFunctionReturn(0);
326fef7b6d8SPeter Brune       }
327fef7b6d8SPeter Brune     }
328e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
329fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3305d10c001SPeter Brune       PetscFunctionReturn(0);
331fef7b6d8SPeter Brune     }
332fef7b6d8SPeter Brune     /* Monitor convergence */
3339566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
334169e2be8SPeter Brune     snes->iter  = i;
335fef7b6d8SPeter Brune     snes->norm  = fnorm;
336c1e67a49SFande Kong     snes->xnorm = xnorm;
337c1e67a49SFande Kong     snes->ynorm = ynorm;
3389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3399566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
3409566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
341302dbbaeSPeter Brune 
342fef7b6d8SPeter Brune     /* Test for convergence */
343dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
3445d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
345302dbbaeSPeter Brune 
346302dbbaeSPeter Brune     /* Call general purpose update function */
347dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
348efd4aadfSBarry Smith     if (snes->npc) {
349a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3509566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, dX));
3519566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
352b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
353b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
354b7281c8aSPeter Brune           PetscFunctionReturn(0);
355b7281c8aSPeter Brune         }
3569566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX, F));
357a71552e2SPeter Brune       } else {
3589566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, dX));
3599566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
360b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
361b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
362b7281c8aSPeter Brune           PetscFunctionReturn(0);
363b7281c8aSPeter Brune         }
364302dbbaeSPeter Brune       }
365a3685310SPeter Brune     } else {
3669566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, dX));
367302dbbaeSPeter Brune     }
368302dbbaeSPeter Brune 
36929c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3700a844d1aSPeter Brune     switch (ncg->type) {
3710a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
37232cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3739566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
37432cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / dXolddotdXold);
375dfb256c7SPeter Brune       break;
3760a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
37732cc618eSPeter Brune       dXolddotdXold = dXdotdX;
3789566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
3799566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
3809566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3819566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
38232cc618eSPeter Brune       beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
383dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
384dfb256c7SPeter Brune       break;
3850a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
3869566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3879566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
3889566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3899566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3909566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3919566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
3929566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3939566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
39432cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
395dfb256c7SPeter Brune       break;
3960a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
3979566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3989566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3999566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4009566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4019566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4029566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4032f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
404dfb256c7SPeter Brune       break;
4050a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
4069566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4079566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4089566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4099566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4102f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
411dfb256c7SPeter Brune       break;
412dfb256c7SPeter Brune     }
413*48a46eb9SPierre Jolivet     if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
4149566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
415fef7b6d8SPeter Brune   }
41663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
417fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
418fef7b6d8SPeter Brune   PetscFunctionReturn(0);
419fef7b6d8SPeter Brune }
420fef7b6d8SPeter Brune 
421fef7b6d8SPeter Brune /*MC
422fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
423fef7b6d8SPeter Brune 
424fef7b6d8SPeter Brune   Level: beginner
425fef7b6d8SPeter Brune 
426fef7b6d8SPeter Brune   Options Database:
4274f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
42841a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
429b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
430fef7b6d8SPeter Brune 
43195452b02SPatrick Sanan    Notes:
43295452b02SPatrick Sanan     This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
433fef7b6d8SPeter Brune           gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
434fef7b6d8SPeter Brune           chooses the initial search direction as F(x) for the initial guess x.
435fef7b6d8SPeter Brune 
4362d547940SBarry Smith           Only supports left non-linear preconditioning.
4372d547940SBarry Smith 
438b5badacbSBarry Smith     Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
439b5badacbSBarry Smith     SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR
440b5badacbSBarry Smith 
4414f02bc6aSBarry Smith    References:
442606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4434f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
444cc3c3c0fSMatthew G. Knepley 
445db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESNCGGetType()`, `SNESLineSearchSetType()`
446fef7b6d8SPeter Brune M*/
4479371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) {
448ea630c6eSPeter Brune   SNES_NCG *neP;
449fef7b6d8SPeter Brune 
450fef7b6d8SPeter Brune   PetscFunctionBegin;
451fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
452fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
453fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
454fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
455fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
456fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
457fef7b6d8SPeter Brune 
458fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
459efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
460efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
461fef7b6d8SPeter Brune 
4624fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4634fc747eaSLawrence Mitchell 
46488976e71SPeter Brune   if (!snes->tolerancesset) {
4650e444f03SPeter Brune     snes->max_funcs = 30000;
4660e444f03SPeter Brune     snes->max_its   = 10000;
467c522fa08SPeter Brune     snes->stol      = 1e-20;
46888976e71SPeter Brune   }
4690e444f03SPeter Brune 
4709566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &neP));
471fef7b6d8SPeter Brune   snes->data   = (void *)neP;
4720298fd71SBarry Smith   neP->monitor = NULL;
4730a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
4749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
475fef7b6d8SPeter Brune   PetscFunctionReturn(0);
476fef7b6d8SPeter Brune }
477