xref: /petsc/src/snes/impls/ncg/snesncg.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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 
4b5badacbSBarry Smith static PetscErrorCode SNESReset_NCG(SNES snes)
5fef7b6d8SPeter Brune {
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 */
18b5badacbSBarry Smith static PetscErrorCode SNESDestroy_NCG(SNES snes)
19fef7b6d8SPeter Brune {
20fef7b6d8SPeter Brune   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
22fef7b6d8SPeter Brune   PetscFunctionReturn(0);
23fef7b6d8SPeter Brune }
24fef7b6d8SPeter Brune 
25fef7b6d8SPeter Brune /*
26fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
27fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
28fef7b6d8SPeter Brune 
29fef7b6d8SPeter Brune    Input Parameters:
30fef7b6d8SPeter Brune +  snes - the SNES context
31fef7b6d8SPeter Brune -  x - the solution vector
32fef7b6d8SPeter Brune 
33fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
34fef7b6d8SPeter Brune  */
35bbd5d0b3SPeter Brune 
36b5badacbSBarry Smith static PetscErrorCode SNESSetUp_NCG(SNES snes)
37fef7b6d8SPeter Brune {
38fef7b6d8SPeter Brune   PetscFunctionBegin;
399566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes,2));
402c71b3e2SJacob Faibussowitsch   PetscCheckFalse(snes->npcside== PC_RIGHT,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
416c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
42fef7b6d8SPeter Brune   PetscFunctionReturn(0);
43fef7b6d8SPeter Brune }
44fef7b6d8SPeter Brune 
45b5badacbSBarry Smith static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
46ea630c6eSPeter Brune {
47ea630c6eSPeter Brune   PetscScalar    alpha, ptAp;
48bbd5d0b3SPeter Brune   Vec            X, Y, F, W;
49bbd5d0b3SPeter Brune   SNES           snes;
50bbd5d0b3SPeter Brune   PetscReal      *fnorm, *xnorm, *ynorm;
51ea630c6eSPeter Brune 
52ea630c6eSPeter Brune   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
54bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
55bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
56bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
57bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
58bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
59bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
60bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
61bbd5d0b3SPeter Brune 
629105210eSPeter Brune   if (!snes->jacobian) {
639566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
649105210eSPeter Brune   }
659105210eSPeter Brune 
66ea630c6eSPeter Brune   /*
67ea630c6eSPeter Brune 
68d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
69ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
70ea630c6eSPeter Brune    */
719566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
729566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &alpha));
739566063dSJacob Faibussowitsch   PetscCall(MatMult(snes->jacobian, Y, W));
749566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, W, &ptAp));
75ea630c6eSPeter Brune   alpha = alpha / ptAp;
769566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,-alpha,Y));
779566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes,X,F));
78bbd5d0b3SPeter Brune 
799566063dSJacob Faibussowitsch   PetscCall(VecNorm(F,NORM_2,fnorm));
809566063dSJacob Faibussowitsch   PetscCall(VecNorm(X,NORM_2,xnorm));
819566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_2,ynorm));
82ea630c6eSPeter Brune   PetscFunctionReturn(0);
83ea630c6eSPeter Brune }
84ea630c6eSPeter Brune 
85b5badacbSBarry Smith /*MC
86b5badacbSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG
87b5badacbSBarry Smith 
88b5badacbSBarry 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.
89b5badacbSBarry 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.
90b5badacbSBarry Smith 
91b5badacbSBarry Smith    Notes: 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 
97b5badacbSBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
98b5badacbSBarry Smith M*/
99b5badacbSBarry Smith 
1008cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
101bbd5d0b3SPeter Brune {
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 */
120b5badacbSBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
121b5badacbSBarry Smith {
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;
1289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(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));
130b5badacbSBarry Smith   if (debug) {
1312f613bf5SBarry Smith     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
132b5badacbSBarry Smith   }
1339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL));
1349566063dSJacob Faibussowitsch   PetscCall(SNESNCGSetType(snes, ncgtype));
1359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
136b5badacbSBarry Smith   if (!snes->linesearch) {
1379566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
138ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
139b5badacbSBarry Smith       if (!snes->npc) {
1409566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
141b5badacbSBarry Smith       } else {
1429566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
143b5badacbSBarry Smith       }
144b5badacbSBarry Smith     }
145ec786807SJed Brown   }
146b5badacbSBarry Smith   PetscFunctionReturn(0);
147b5badacbSBarry Smith }
148b5badacbSBarry Smith 
149b5badacbSBarry Smith /*
150b5badacbSBarry Smith   SNESView_NCG - Prints info from the SNESNCG data structure.
151b5badacbSBarry Smith 
152b5badacbSBarry Smith   Input Parameters:
153b5badacbSBarry Smith + SNES - the SNES context
154b5badacbSBarry Smith - viewer - visualization context
155b5badacbSBarry Smith 
156b5badacbSBarry Smith   Application Interface Routine: SNESView()
157b5badacbSBarry Smith */
158b5badacbSBarry Smith static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
159b5badacbSBarry Smith {
160b5badacbSBarry Smith   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
161b5badacbSBarry Smith   PetscBool      iascii;
162b5badacbSBarry Smith 
163b5badacbSBarry Smith   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
165b5badacbSBarry Smith   if (iascii) {
1669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
167b5badacbSBarry Smith   }
168b5badacbSBarry Smith   PetscFunctionReturn(0);
169b5badacbSBarry Smith }
170b5badacbSBarry Smith 
171b5badacbSBarry Smith /*
17288d374b2SPeter Brune 
17388d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
17488d374b2SPeter Brune 
175b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
176b5badacbSBarry Smith 
177b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
178b5badacbSBarry Smith 
17988d374b2SPeter Brune  */
18067392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
18167392de3SBarry Smith {
18288d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
18367392de3SBarry Smith 
18488d374b2SPeter Brune   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &ftf));
1869566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty));
18788d374b2SPeter Brune   h      = 1e-5*fty / fty;
1889566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
1899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W, -h, Y));          /* this is arbitrary */
1909566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, W, G));
1919566063dSJacob Faibussowitsch   PetscCall(VecDot(G, F, &ftg));
19288d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
19388d374b2SPeter Brune   PetscFunctionReturn(0);
19488d374b2SPeter Brune }
19588d374b2SPeter Brune 
1960a844d1aSPeter Brune /*@
1970a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
1980a844d1aSPeter Brune 
1990a844d1aSPeter Brune     Logically Collective on SNES
2000a844d1aSPeter Brune 
2010a844d1aSPeter Brune     Input Parameters:
2020a844d1aSPeter Brune +   snes - the iterative context
2030a844d1aSPeter Brune -   btype - update type
2040a844d1aSPeter Brune 
2050a844d1aSPeter Brune     Options Database:
206b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
2070a844d1aSPeter Brune 
2080a844d1aSPeter Brune     Level: intermediate
2090a844d1aSPeter Brune 
2100a844d1aSPeter Brune     SNESNCGSelectTypes:
2110a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2120a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2130a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2140a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2150a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2160a844d1aSPeter Brune 
2170a844d1aSPeter Brune    Notes:
218b5badacbSBarry Smith    SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.
219b5badacbSBarry Smith 
220b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
221b5badacbSBarry Smith    that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?
2220a844d1aSPeter Brune 
2230a844d1aSPeter Brune @*/
2240adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
2250adebc6cSBarry Smith {
2260a844d1aSPeter Brune   PetscFunctionBegin;
2270a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
228*cac4c232SBarry Smith   PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
2290a844d1aSPeter Brune   PetscFunctionReturn(0);
2300a844d1aSPeter Brune }
2310a844d1aSPeter Brune 
232b5badacbSBarry Smith static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
2330adebc6cSBarry Smith {
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 */
252b5badacbSBarry Smith static PetscErrorCode SNESSolve_NCG(SNES snes)
253fef7b6d8SPeter Brune {
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;
2642c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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 */
3249566063dSJacob Faibussowitsch   PetscCall((*snes->ops->converged)(snes,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 */
328fef7b6d8SPeter Brune   if (snes->ops->update) {
3299566063dSJacob Faibussowitsch     PetscCall((*snes->ops->update)(snes, snes->iter));
330fef7b6d8SPeter Brune   }
331fef7b6d8SPeter Brune 
332fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
333bbd5d0b3SPeter Brune 
33409c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
33529c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3360a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
3379566063dSJacob Faibussowitsch       PetscCall(VecCopy(dX, dXold));
33888d374b2SPeter Brune     }
3399566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch,X,F,&fnorm,lX));
3409566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
3419566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
342422a814eSBarry Smith     if (lsresult) {
343fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
344fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3455d10c001SPeter Brune         PetscFunctionReturn(0);
346fef7b6d8SPeter Brune       }
347fef7b6d8SPeter Brune     }
348e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
349fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3505d10c001SPeter Brune       PetscFunctionReturn(0);
351fef7b6d8SPeter Brune     }
352fef7b6d8SPeter Brune     /* Monitor convergence */
3539566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
354169e2be8SPeter Brune     snes->iter = i;
355fef7b6d8SPeter Brune     snes->norm = fnorm;
356c1e67a49SFande Kong     snes->xnorm = xnorm;
357c1e67a49SFande Kong     snes->ynorm = ynorm;
3589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3599566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
3609566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
361302dbbaeSPeter Brune 
362fef7b6d8SPeter Brune     /* Test for convergence */
3639566063dSJacob Faibussowitsch     PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP));
3645d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
365302dbbaeSPeter Brune 
366302dbbaeSPeter Brune     /* Call general purpose update function */
367302dbbaeSPeter Brune     if (snes->ops->update) {
3689566063dSJacob Faibussowitsch       PetscCall((*snes->ops->update)(snes, snes->iter));
369302dbbaeSPeter Brune     }
370efd4aadfSBarry Smith     if (snes->npc) {
371a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3729566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes,X,NULL,dX));
3739566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc,&reason));
374b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
375b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
376b7281c8aSPeter Brune           PetscFunctionReturn(0);
377b7281c8aSPeter Brune         }
3789566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX,F));
379a71552e2SPeter Brune       } else {
3809566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes,X,F,dX));
3819566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc,&reason));
382b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
383b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
384b7281c8aSPeter Brune           PetscFunctionReturn(0);
385b7281c8aSPeter Brune         }
386302dbbaeSPeter Brune       }
387a3685310SPeter Brune     } else {
3889566063dSJacob Faibussowitsch       PetscCall(VecCopy(F,dX));
389302dbbaeSPeter Brune     }
390302dbbaeSPeter Brune 
39129c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3920a844d1aSPeter Brune     switch (ncg->type) {
3930a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
39432cc618eSPeter Brune       dXolddotdXold= dXdotdX;
3959566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
39632cc618eSPeter Brune       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
397dfb256c7SPeter Brune       break;
3980a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
39932cc618eSPeter Brune       dXolddotdXold= dXdotdX;
4009566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
4019566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
4029566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4039566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
40432cc618eSPeter Brune       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
405dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
406dfb256c7SPeter Brune       break;
4070a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
4089566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4099566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
4109566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4119566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4129566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4139566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
4149566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4159566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
41632cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
417dfb256c7SPeter Brune       break;
4180a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
4199566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4209566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4219566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4229566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4239566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4249566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4252f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
426dfb256c7SPeter Brune       break;
4270a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
4289566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4299566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4309566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4319566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4322f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
433dfb256c7SPeter Brune       break;
434dfb256c7SPeter Brune     }
435dfb256c7SPeter Brune     if (ncg->monitor) {
4369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
437dfb256c7SPeter Brune     }
4389566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
439fef7b6d8SPeter Brune   }
4409566063dSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %D\n", maxits));
441fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
442fef7b6d8SPeter Brune   PetscFunctionReturn(0);
443fef7b6d8SPeter Brune }
444fef7b6d8SPeter Brune 
445fef7b6d8SPeter Brune /*MC
446fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
447fef7b6d8SPeter Brune 
448fef7b6d8SPeter Brune   Level: beginner
449fef7b6d8SPeter Brune 
450fef7b6d8SPeter Brune   Options Database:
4514f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
45241a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
453b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
454fef7b6d8SPeter Brune 
45595452b02SPatrick Sanan    Notes:
45695452b02SPatrick Sanan     This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
457fef7b6d8SPeter Brune           gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
458fef7b6d8SPeter Brune           chooses the initial search direction as F(x) for the initial guess x.
459fef7b6d8SPeter Brune 
4602d547940SBarry Smith           Only supports left non-linear preconditioning.
4612d547940SBarry Smith 
462b5badacbSBarry Smith     Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
463b5badacbSBarry Smith     SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR
464b5badacbSBarry Smith 
4654f02bc6aSBarry Smith    References:
466606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4674f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
468cc3c3c0fSMatthew G. Knepley 
469b5badacbSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType()
470fef7b6d8SPeter Brune M*/
4718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
472fef7b6d8SPeter Brune {
473ea630c6eSPeter Brune   SNES_NCG       * neP;
474fef7b6d8SPeter Brune 
475fef7b6d8SPeter Brune   PetscFunctionBegin;
476fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
477fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
478fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
479fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
480fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
481fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
482fef7b6d8SPeter Brune 
483fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
484efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
485efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
486fef7b6d8SPeter Brune 
4874fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4884fc747eaSLawrence Mitchell 
48988976e71SPeter Brune   if (!snes->tolerancesset) {
4900e444f03SPeter Brune     snes->max_funcs = 30000;
4910e444f03SPeter Brune     snes->max_its   = 10000;
492c522fa08SPeter Brune     snes->stol      = 1e-20;
49388976e71SPeter Brune   }
4940e444f03SPeter Brune 
4959566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&neP));
496fef7b6d8SPeter Brune   snes->data   = (void*) neP;
4970298fd71SBarry Smith   neP->monitor = NULL;
4980a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG));
500fef7b6d8SPeter Brune   PetscFunctionReturn(0);
501fef7b6d8SPeter Brune }
502