xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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;
21*2e956fe4SStefano 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 
37b5badacbSBarry Smith static PetscErrorCode SNESSetUp_NCG(SNES snes)
38fef7b6d8SPeter Brune {
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 
46b5badacbSBarry Smith static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
47ea630c6eSPeter Brune {
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 
639105210eSPeter Brune   if (!snes->jacobian) {
649566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
659105210eSPeter Brune   }
669105210eSPeter Brune 
67ea630c6eSPeter Brune   /*
68ea630c6eSPeter Brune 
69d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
70ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
71ea630c6eSPeter Brune    */
729566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
739566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &alpha));
749566063dSJacob Faibussowitsch   PetscCall(MatMult(snes->jacobian, Y, W));
759566063dSJacob Faibussowitsch   PetscCall(VecDot(Y, W, &ptAp));
76ea630c6eSPeter Brune   alpha = alpha / ptAp;
779566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,-alpha,Y));
789566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes,X,F));
79bbd5d0b3SPeter Brune 
809566063dSJacob Faibussowitsch   PetscCall(VecNorm(F,NORM_2,fnorm));
819566063dSJacob Faibussowitsch   PetscCall(VecNorm(X,NORM_2,xnorm));
829566063dSJacob Faibussowitsch   PetscCall(VecNorm(Y,NORM_2,ynorm));
83ea630c6eSPeter Brune   PetscFunctionReturn(0);
84ea630c6eSPeter Brune }
85ea630c6eSPeter Brune 
86b5badacbSBarry Smith /*MC
87b5badacbSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG
88b5badacbSBarry Smith 
89b5badacbSBarry 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.
90b5badacbSBarry 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.
91b5badacbSBarry Smith 
92b5badacbSBarry Smith    Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
93b5badacbSBarry Smith 
94b5badacbSBarry Smith    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
95b5badacbSBarry Smith 
96b5badacbSBarry Smith    Level: advanced
97b5badacbSBarry Smith 
98db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
99b5badacbSBarry Smith M*/
100b5badacbSBarry Smith 
1018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
102bbd5d0b3SPeter Brune {
103bbd5d0b3SPeter Brune   PetscFunctionBegin;
104f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
1050298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1060298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1070298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1080298fd71SBarry Smith   linesearch->ops->view           = NULL;
1090298fd71SBarry Smith   linesearch->ops->setup          = NULL;
110bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
111bbd5d0b3SPeter Brune }
112bbd5d0b3SPeter Brune 
11388d374b2SPeter Brune /*
114b5badacbSBarry Smith   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
115b5badacbSBarry Smith 
116b5badacbSBarry Smith   Input Parameter:
117b5badacbSBarry Smith . snes - the SNES context
118b5badacbSBarry Smith 
119b5badacbSBarry Smith   Application Interface Routine: SNESSetFromOptions()
120b5badacbSBarry Smith */
121b5badacbSBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
122b5badacbSBarry Smith {
123b5badacbSBarry Smith   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
124b5badacbSBarry Smith   PetscBool      debug = PETSC_FALSE;
125b5badacbSBarry Smith   SNESNCGType    ncgtype=ncg->type;
126b5badacbSBarry Smith   SNESLineSearch linesearch;
127b5badacbSBarry Smith 
128b5badacbSBarry Smith   PetscFunctionBegin;
129d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SNES NCG options");
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
131b5badacbSBarry Smith   if (debug) {
1322f613bf5SBarry Smith     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
133b5badacbSBarry Smith   }
1349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL));
1359566063dSJacob Faibussowitsch   PetscCall(SNESNCGSetType(snes, ncgtype));
136d0609cedSBarry Smith   PetscOptionsHeadEnd();
137b5badacbSBarry Smith   if (!snes->linesearch) {
1389566063dSJacob Faibussowitsch     PetscCall(SNESGetLineSearch(snes, &linesearch));
139ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
140b5badacbSBarry Smith       if (!snes->npc) {
1419566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
142b5badacbSBarry Smith       } else {
1439566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
144b5badacbSBarry Smith       }
145b5badacbSBarry Smith     }
146ec786807SJed Brown   }
147b5badacbSBarry Smith   PetscFunctionReturn(0);
148b5badacbSBarry Smith }
149b5badacbSBarry Smith 
150b5badacbSBarry Smith /*
151b5badacbSBarry Smith   SNESView_NCG - Prints info from the SNESNCG data structure.
152b5badacbSBarry Smith 
153b5badacbSBarry Smith   Input Parameters:
154b5badacbSBarry Smith + SNES - the SNES context
155b5badacbSBarry Smith - viewer - visualization context
156b5badacbSBarry Smith 
157b5badacbSBarry Smith   Application Interface Routine: SNESView()
158b5badacbSBarry Smith */
159b5badacbSBarry Smith static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
160b5badacbSBarry Smith {
161b5badacbSBarry Smith   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
162b5badacbSBarry Smith   PetscBool      iascii;
163b5badacbSBarry Smith 
164b5badacbSBarry Smith   PetscFunctionBegin;
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
166b5badacbSBarry Smith   if (iascii) {
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
168b5badacbSBarry Smith   }
169b5badacbSBarry Smith   PetscFunctionReturn(0);
170b5badacbSBarry Smith }
171b5badacbSBarry Smith 
172b5badacbSBarry Smith /*
17388d374b2SPeter Brune 
17488d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
17588d374b2SPeter Brune 
176b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
177b5badacbSBarry Smith 
178b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
179b5badacbSBarry Smith 
18088d374b2SPeter Brune  */
18167392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
18267392de3SBarry Smith {
18388d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
18467392de3SBarry Smith 
18588d374b2SPeter Brune   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(VecDot(F, F, &ftf));
1879566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty));
18888d374b2SPeter Brune   h      = 1e-5*fty / fty;
1899566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, W));
1909566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W, -h, Y));          /* this is arbitrary */
1919566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, W, G));
1929566063dSJacob Faibussowitsch   PetscCall(VecDot(G, F, &ftg));
19388d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
19488d374b2SPeter Brune   PetscFunctionReturn(0);
19588d374b2SPeter Brune }
19688d374b2SPeter Brune 
1970a844d1aSPeter Brune /*@
1980a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
1990a844d1aSPeter Brune 
2000a844d1aSPeter Brune     Logically Collective on SNES
2010a844d1aSPeter Brune 
2020a844d1aSPeter Brune     Input Parameters:
2030a844d1aSPeter Brune +   snes - the iterative context
2040a844d1aSPeter Brune -   btype - update type
2050a844d1aSPeter Brune 
2060a844d1aSPeter Brune     Options Database:
207b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
2080a844d1aSPeter Brune 
2090a844d1aSPeter Brune     Level: intermediate
2100a844d1aSPeter Brune 
2110a844d1aSPeter Brune     SNESNCGSelectTypes:
2120a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2130a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2140a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2150a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2160a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2170a844d1aSPeter Brune 
2180a844d1aSPeter Brune    Notes:
219b5badacbSBarry Smith    SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.
220b5badacbSBarry Smith 
221b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
222b5badacbSBarry Smith    that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?
2230a844d1aSPeter Brune 
2240a844d1aSPeter Brune @*/
2250adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
2260adebc6cSBarry Smith {
2270a844d1aSPeter Brune   PetscFunctionBegin;
2280a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
229cac4c232SBarry Smith   PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
2300a844d1aSPeter Brune   PetscFunctionReturn(0);
2310a844d1aSPeter Brune }
2320a844d1aSPeter Brune 
233b5badacbSBarry Smith static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
2340adebc6cSBarry Smith {
2350a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG*)snes->data;
2366e111a19SKarl Rupp 
2370a844d1aSPeter Brune   PetscFunctionBegin;
2380a844d1aSPeter Brune   ncg->type = btype;
2390a844d1aSPeter Brune   PetscFunctionReturn(0);
2400a844d1aSPeter Brune }
2410a844d1aSPeter Brune 
242fef7b6d8SPeter Brune /*
243fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
244fef7b6d8SPeter Brune 
245fef7b6d8SPeter Brune   Input Parameters:
246fef7b6d8SPeter Brune . snes - the SNES context
247fef7b6d8SPeter Brune 
248fef7b6d8SPeter Brune   Output Parameter:
249fef7b6d8SPeter Brune . outits - number of iterations until termination
250fef7b6d8SPeter Brune 
251fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
252fef7b6d8SPeter Brune */
253b5badacbSBarry Smith static PetscErrorCode SNESSolve_NCG(SNES snes)
254fef7b6d8SPeter Brune {
255dfb256c7SPeter Brune   SNES_NCG             *ncg = (SNES_NCG*)snes->data;
25632cc618eSPeter Brune   Vec                  X,dX,lX,F,dXold;
257bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
25832cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
259fef7b6d8SPeter Brune   PetscInt             maxits, i;
260422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
261f1c6b773SPeter Brune   SNESLineSearch       linesearch;
262b7281c8aSPeter Brune   SNESConvergedReason  reason;
26388d374b2SPeter Brune 
264fef7b6d8SPeter Brune   PetscFunctionBegin;
2650b121fc5SBarry 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);
266c579b300SPatrick Farrell 
2679566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite));
268fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
269fef7b6d8SPeter Brune 
270fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
271fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
27232cc618eSPeter Brune   dXold  = snes->work[0];            /* The previous iterate of X */
27388d374b2SPeter Brune   dX     = snes->work[1];            /* the preconditioned direction */
274169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
275fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
276fef7b6d8SPeter Brune 
2779566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
2789e764e56SPeter Brune 
2799566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
280fef7b6d8SPeter Brune   snes->iter = 0;
281fef7b6d8SPeter Brune   snes->norm = 0.;
2829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
283fef7b6d8SPeter Brune 
284bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
285a71552e2SPeter Brune 
286efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2879566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes,X,NULL,dX));
2889566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc,&reason));
289b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
290b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
291b7281c8aSPeter Brune       PetscFunctionReturn(0);
292b7281c8aSPeter Brune     }
2939566063dSJacob Faibussowitsch     PetscCall(VecCopy(dX,F));
2949566063dSJacob Faibussowitsch     PetscCall(VecNorm(F,NORM_2,&fnorm));
295a71552e2SPeter Brune   } else {
296e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2979566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes,X,F));
2981aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
299c1c75074SPeter Brune 
300fef7b6d8SPeter Brune     /* convergence test */
3019566063dSJacob Faibussowitsch     PetscCall(VecNorm(F,NORM_2,&fnorm));
302422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
3039566063dSJacob Faibussowitsch     PetscCall(VecCopy(F,dX));
304a71552e2SPeter Brune   }
305efd4aadfSBarry Smith   if (snes->npc) {
306a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
3079566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes,X,F,dX));
3089566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc,&reason));
309b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
310b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
311b7281c8aSPeter Brune         PetscFunctionReturn(0);
312b7281c8aSPeter Brune       }
313a71552e2SPeter Brune     }
314a71552e2SPeter Brune   }
3159566063dSJacob Faibussowitsch   PetscCall(VecCopy(dX,lX));
3169566063dSJacob Faibussowitsch   PetscCall(VecDot(dX, dX, &dXdotdX));
317a71552e2SPeter Brune 
3189566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
319fef7b6d8SPeter Brune   snes->norm = fnorm;
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3219566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes,fnorm,0));
3229566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes,0,fnorm));
323fef7b6d8SPeter Brune 
324fef7b6d8SPeter Brune   /* test convergence */
3259566063dSJacob Faibussowitsch   PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP));
326fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
327fef7b6d8SPeter Brune 
328fef7b6d8SPeter Brune   /* Call general purpose update function */
329fef7b6d8SPeter Brune   if (snes->ops->update) {
3309566063dSJacob Faibussowitsch     PetscCall((*snes->ops->update)(snes, snes->iter));
331fef7b6d8SPeter Brune   }
332fef7b6d8SPeter Brune 
333fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
334bbd5d0b3SPeter Brune 
33509c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
33629c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3370a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
3389566063dSJacob Faibussowitsch       PetscCall(VecCopy(dX, dXold));
33988d374b2SPeter Brune     }
3409566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch,X,F,&fnorm,lX));
3419566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
3429566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
343422a814eSBarry Smith     if (lsresult) {
344fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
345fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3465d10c001SPeter Brune         PetscFunctionReturn(0);
347fef7b6d8SPeter Brune       }
348fef7b6d8SPeter Brune     }
349e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
350fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3515d10c001SPeter Brune       PetscFunctionReturn(0);
352fef7b6d8SPeter Brune     }
353fef7b6d8SPeter Brune     /* Monitor convergence */
3549566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
355169e2be8SPeter Brune     snes->iter = i;
356fef7b6d8SPeter Brune     snes->norm = fnorm;
357c1e67a49SFande Kong     snes->xnorm = xnorm;
358c1e67a49SFande Kong     snes->ynorm = ynorm;
3599566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3609566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
3619566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
362302dbbaeSPeter Brune 
363fef7b6d8SPeter Brune     /* Test for convergence */
3649566063dSJacob Faibussowitsch     PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP));
3655d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
366302dbbaeSPeter Brune 
367302dbbaeSPeter Brune     /* Call general purpose update function */
368302dbbaeSPeter Brune     if (snes->ops->update) {
3699566063dSJacob Faibussowitsch       PetscCall((*snes->ops->update)(snes, snes->iter));
370302dbbaeSPeter Brune     }
371efd4aadfSBarry Smith     if (snes->npc) {
372a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3739566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes,X,NULL,dX));
3749566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc,&reason));
375b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
376b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
377b7281c8aSPeter Brune           PetscFunctionReturn(0);
378b7281c8aSPeter Brune         }
3799566063dSJacob Faibussowitsch         PetscCall(VecCopy(dX,F));
380a71552e2SPeter Brune       } else {
3819566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes,X,F,dX));
3829566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc,&reason));
383b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
384b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
385b7281c8aSPeter Brune           PetscFunctionReturn(0);
386b7281c8aSPeter Brune         }
387302dbbaeSPeter Brune       }
388a3685310SPeter Brune     } else {
3899566063dSJacob Faibussowitsch       PetscCall(VecCopy(F,dX));
390302dbbaeSPeter Brune     }
391302dbbaeSPeter Brune 
39229c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3930a844d1aSPeter Brune     switch (ncg->type) {
3940a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
39532cc618eSPeter Brune       dXolddotdXold= dXdotdX;
3969566063dSJacob Faibussowitsch       PetscCall(VecDot(dX, dX, &dXdotdX));
39732cc618eSPeter Brune       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
398dfb256c7SPeter Brune       break;
3990a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
40032cc618eSPeter Brune       dXolddotdXold= dXdotdX;
4019566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdXold));
4029566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
4039566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4049566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
40532cc618eSPeter Brune       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
406dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
407dfb256c7SPeter Brune       break;
4080a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
4099566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4109566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
4119566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4129566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4139566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4149566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
4159566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4169566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
41732cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
418dfb256c7SPeter Brune       break;
4190a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
4209566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4219566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
4229566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4239566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4249566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
4259566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4262f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
427dfb256c7SPeter Brune       break;
4280a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
4299566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
4309566063dSJacob Faibussowitsch       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
4319566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
4329566063dSJacob Faibussowitsch       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
4332f613bf5SBarry Smith       beta = PetscRealPart(dXdotdX / lXdotdXold);
434dfb256c7SPeter Brune       break;
435dfb256c7SPeter Brune     }
436dfb256c7SPeter Brune     if (ncg->monitor) {
4379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
438dfb256c7SPeter Brune     }
4399566063dSJacob Faibussowitsch     PetscCall(VecAYPX(lX, beta, dX));
440fef7b6d8SPeter Brune   }
44163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
442fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
443fef7b6d8SPeter Brune   PetscFunctionReturn(0);
444fef7b6d8SPeter Brune }
445fef7b6d8SPeter Brune 
446fef7b6d8SPeter Brune /*MC
447fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
448fef7b6d8SPeter Brune 
449fef7b6d8SPeter Brune   Level: beginner
450fef7b6d8SPeter Brune 
451fef7b6d8SPeter Brune   Options Database:
4524f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
45341a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
454b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
455fef7b6d8SPeter Brune 
45695452b02SPatrick Sanan    Notes:
45795452b02SPatrick Sanan     This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
458fef7b6d8SPeter Brune           gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
459fef7b6d8SPeter Brune           chooses the initial search direction as F(x) for the initial guess x.
460fef7b6d8SPeter Brune 
4612d547940SBarry Smith           Only supports left non-linear preconditioning.
4622d547940SBarry Smith 
463b5badacbSBarry Smith     Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
464b5badacbSBarry Smith     SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR
465b5badacbSBarry Smith 
4664f02bc6aSBarry Smith    References:
467606c0280SSatish Balay .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4684f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
469cc3c3c0fSMatthew G. Knepley 
470db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESNCGGetType()`, `SNESLineSearchSetType()`
471fef7b6d8SPeter Brune M*/
4728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
473fef7b6d8SPeter Brune {
474ea630c6eSPeter Brune   SNES_NCG       * neP;
475fef7b6d8SPeter Brune 
476fef7b6d8SPeter Brune   PetscFunctionBegin;
477fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
478fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
479fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
480fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
481fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
482fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
483fef7b6d8SPeter Brune 
484fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
485efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
486efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
487fef7b6d8SPeter Brune 
4884fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
4894fc747eaSLawrence Mitchell 
49088976e71SPeter Brune   if (!snes->tolerancesset) {
4910e444f03SPeter Brune     snes->max_funcs = 30000;
4920e444f03SPeter Brune     snes->max_its   = 10000;
493c522fa08SPeter Brune     snes->stol      = 1e-20;
49488976e71SPeter Brune   }
4950e444f03SPeter Brune 
4969566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&neP));
497fef7b6d8SPeter Brune   snes->data   = (void*) neP;
4980298fd71SBarry Smith   neP->monitor = NULL;
4990a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG));
501fef7b6d8SPeter Brune   PetscFunctionReturn(0);
502fef7b6d8SPeter Brune }
503