xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 9e5d08929599755f0fc0004c3575ee121dac326a)
10a844d1aSPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/
2*9e5d0892SLisandro 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   PetscErrorCode ierr;
21fef7b6d8SPeter Brune 
22fef7b6d8SPeter Brune   PetscFunctionBegin;
23fef7b6d8SPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
24fef7b6d8SPeter Brune   PetscFunctionReturn(0);
25fef7b6d8SPeter Brune }
26fef7b6d8SPeter Brune 
27fef7b6d8SPeter Brune /*
28fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
29fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
30fef7b6d8SPeter Brune 
31fef7b6d8SPeter Brune    Input Parameters:
32fef7b6d8SPeter Brune +  snes - the SNES context
33fef7b6d8SPeter Brune -  x - the solution vector
34fef7b6d8SPeter Brune 
35fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
36fef7b6d8SPeter Brune  */
37bbd5d0b3SPeter Brune 
38b5badacbSBarry Smith static PetscErrorCode SNESSetUp_NCG(SNES snes)
39fef7b6d8SPeter Brune {
40fef7b6d8SPeter Brune   PetscErrorCode ierr;
41fef7b6d8SPeter Brune 
42fef7b6d8SPeter Brune   PetscFunctionBegin;
43fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr);
44efd4aadfSBarry Smith   if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
456c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
46fef7b6d8SPeter Brune   PetscFunctionReturn(0);
47fef7b6d8SPeter Brune }
48fef7b6d8SPeter Brune 
49b5badacbSBarry Smith static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
50ea630c6eSPeter Brune {
51ea630c6eSPeter Brune   PetscScalar    alpha, ptAp;
52bbd5d0b3SPeter Brune   Vec            X, Y, F, W;
53bbd5d0b3SPeter Brune   SNES           snes;
54ea630c6eSPeter Brune   PetscErrorCode ierr;
55bbd5d0b3SPeter Brune   PetscReal      *fnorm, *xnorm, *ynorm;
56ea630c6eSPeter Brune 
57ea630c6eSPeter Brune   PetscFunctionBegin;
58f1c6b773SPeter Brune   ierr  = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
59bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
60bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
61bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
62bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
63bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
64bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
65bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
66bbd5d0b3SPeter Brune 
679105210eSPeter Brune   if (!snes->jacobian) {
689105210eSPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
699105210eSPeter Brune   }
709105210eSPeter Brune 
71ea630c6eSPeter Brune   /*
72ea630c6eSPeter Brune 
73d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
74ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
75ea630c6eSPeter Brune    */
76d1e9a80fSBarry Smith   ierr  = SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);CHKERRQ(ierr);
77ea630c6eSPeter Brune   ierr  = VecDot(F, F, &alpha);CHKERRQ(ierr);
78ea630c6eSPeter Brune   ierr  = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
79ea630c6eSPeter Brune   ierr  = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
80ea630c6eSPeter Brune   alpha = alpha / ptAp;
819105210eSPeter Brune   ierr  = VecAXPY(X,-alpha,Y);CHKERRQ(ierr);
82bbd5d0b3SPeter Brune   ierr  = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
83bbd5d0b3SPeter Brune 
84bbd5d0b3SPeter Brune   ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
85bbd5d0b3SPeter Brune   ierr = VecNorm(X,NORM_2,xnorm);CHKERRQ(ierr);
86bbd5d0b3SPeter Brune   ierr = VecNorm(Y,NORM_2,ynorm);CHKERRQ(ierr);
87ea630c6eSPeter Brune   PetscFunctionReturn(0);
88ea630c6eSPeter Brune }
89ea630c6eSPeter Brune 
90b5badacbSBarry Smith /*MC
91b5badacbSBarry Smith    SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG
92b5badacbSBarry Smith 
93b5badacbSBarry 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.
94b5badacbSBarry 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.
95b5badacbSBarry Smith 
96b5badacbSBarry Smith    Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
97b5badacbSBarry Smith 
98b5badacbSBarry Smith    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
99b5badacbSBarry Smith 
100b5badacbSBarry Smith    Level: advanced
101b5badacbSBarry Smith 
102b5badacbSBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
103b5badacbSBarry Smith M*/
104b5badacbSBarry Smith 
1058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
106bbd5d0b3SPeter Brune {
107bbd5d0b3SPeter Brune   PetscFunctionBegin;
108f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
1090298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1100298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1110298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1120298fd71SBarry Smith   linesearch->ops->view           = NULL;
1130298fd71SBarry Smith   linesearch->ops->setup          = NULL;
114bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
115bbd5d0b3SPeter Brune }
116bbd5d0b3SPeter Brune 
11788d374b2SPeter Brune /*
118b5badacbSBarry Smith   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
119b5badacbSBarry Smith 
120b5badacbSBarry Smith   Input Parameter:
121b5badacbSBarry Smith . snes - the SNES context
122b5badacbSBarry Smith 
123b5badacbSBarry Smith   Application Interface Routine: SNESSetFromOptions()
124b5badacbSBarry Smith */
125b5badacbSBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
126b5badacbSBarry Smith {
127b5badacbSBarry Smith   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
128b5badacbSBarry Smith   PetscErrorCode ierr;
129b5badacbSBarry Smith   PetscBool      debug = PETSC_FALSE;
130b5badacbSBarry Smith   SNESNCGType    ncgtype=ncg->type;
131b5badacbSBarry Smith   SNESLineSearch linesearch;
132b5badacbSBarry Smith 
133b5badacbSBarry Smith   PetscFunctionBegin;
134b5badacbSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES NCG options");CHKERRQ(ierr);
135b5badacbSBarry Smith   ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr);
136b5badacbSBarry Smith   if (debug) {
137b5badacbSBarry Smith     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
138b5badacbSBarry Smith   }
139b5badacbSBarry Smith   ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr);
140b5badacbSBarry Smith   ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr);
141b5badacbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
142b5badacbSBarry Smith   if (!snes->linesearch) {
143b5badacbSBarry Smith     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
144ec786807SJed Brown     if (!((PetscObject)linesearch)->type_name) {
145b5badacbSBarry Smith       if (!snes->npc) {
146b5badacbSBarry Smith         ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
147b5badacbSBarry Smith       } else {
148b5badacbSBarry Smith         ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
149b5badacbSBarry Smith       }
150b5badacbSBarry Smith     }
151ec786807SJed Brown   }
152b5badacbSBarry Smith   PetscFunctionReturn(0);
153b5badacbSBarry Smith }
154b5badacbSBarry Smith 
155b5badacbSBarry Smith /*
156b5badacbSBarry Smith   SNESView_NCG - Prints info from the SNESNCG data structure.
157b5badacbSBarry Smith 
158b5badacbSBarry Smith   Input Parameters:
159b5badacbSBarry Smith + SNES - the SNES context
160b5badacbSBarry Smith - viewer - visualization context
161b5badacbSBarry Smith 
162b5badacbSBarry Smith   Application Interface Routine: SNESView()
163b5badacbSBarry Smith */
164b5badacbSBarry Smith static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
165b5badacbSBarry Smith {
166b5badacbSBarry Smith   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
167b5badacbSBarry Smith   PetscBool      iascii;
168b5badacbSBarry Smith   PetscErrorCode ierr;
169b5badacbSBarry Smith 
170b5badacbSBarry Smith   PetscFunctionBegin;
171b5badacbSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
172b5badacbSBarry Smith   if (iascii) {
173b5badacbSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]);CHKERRQ(ierr);
174b5badacbSBarry Smith   }
175b5badacbSBarry Smith   PetscFunctionReturn(0);
176b5badacbSBarry Smith }
177b5badacbSBarry Smith 
178b5badacbSBarry Smith /*
17988d374b2SPeter Brune 
18088d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
18188d374b2SPeter Brune 
182b5badacbSBarry Smith  This routine is not currently used. I don't know what its intended purpose is.
183b5badacbSBarry Smith 
184b5badacbSBarry Smith  Note it has a hardwired differencing parameter of 1e-5
185b5badacbSBarry Smith 
18688d374b2SPeter Brune  */
18767392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
18867392de3SBarry Smith {
18988d374b2SPeter Brune   PetscErrorCode ierr;
19088d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
19167392de3SBarry Smith 
19288d374b2SPeter Brune   PetscFunctionBegin;
19388d374b2SPeter Brune   ierr   = VecDot(F, F, &ftf);CHKERRQ(ierr);
19488d374b2SPeter Brune   ierr   = VecDot(F, Y, &fty);CHKERRQ(ierr);
19588d374b2SPeter Brune   h      = 1e-5*fty / fty;
19688d374b2SPeter Brune   ierr   = VecCopy(X, W);CHKERRQ(ierr);
19788d374b2SPeter Brune   ierr   = VecAXPY(W, -h, Y);CHKERRQ(ierr);          /* this is arbitrary */
19888d374b2SPeter Brune   ierr   = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
19988d374b2SPeter Brune   ierr   = VecDot(G, F, &ftg);CHKERRQ(ierr);
20088d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
20188d374b2SPeter Brune   PetscFunctionReturn(0);
20288d374b2SPeter Brune }
20388d374b2SPeter Brune 
2040a844d1aSPeter Brune /*@
2050a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
2060a844d1aSPeter Brune 
2070a844d1aSPeter Brune     Logically Collective on SNES
2080a844d1aSPeter Brune 
2090a844d1aSPeter Brune     Input Parameters:
2100a844d1aSPeter Brune +   snes - the iterative context
2110a844d1aSPeter Brune -   btype - update type
2120a844d1aSPeter Brune 
2130a844d1aSPeter Brune     Options Database:
214b5badacbSBarry Smith .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta
2150a844d1aSPeter Brune 
2160a844d1aSPeter Brune     Level: intermediate
2170a844d1aSPeter Brune 
2180a844d1aSPeter Brune     SNESNCGSelectTypes:
2190a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2200a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2210a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2220a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2230a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2240a844d1aSPeter Brune 
2250a844d1aSPeter Brune    Notes:
226b5badacbSBarry Smith    SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.
227b5badacbSBarry Smith 
228b5badacbSBarry Smith    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
229b5badacbSBarry Smith    that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?
2300a844d1aSPeter Brune 
2310a844d1aSPeter Brune @*/
2320adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
2330adebc6cSBarry Smith {
2340a844d1aSPeter Brune   PetscErrorCode ierr;
2356e111a19SKarl Rupp 
2360a844d1aSPeter Brune   PetscFunctionBegin;
2370a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2380a844d1aSPeter Brune   ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr);
2390a844d1aSPeter Brune   PetscFunctionReturn(0);
2400a844d1aSPeter Brune }
2410a844d1aSPeter Brune 
242b5badacbSBarry Smith static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
2430adebc6cSBarry Smith {
2440a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG*)snes->data;
2456e111a19SKarl Rupp 
2460a844d1aSPeter Brune   PetscFunctionBegin;
2470a844d1aSPeter Brune   ncg->type = btype;
2480a844d1aSPeter Brune   PetscFunctionReturn(0);
2490a844d1aSPeter Brune }
2500a844d1aSPeter Brune 
251fef7b6d8SPeter Brune /*
252fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
253fef7b6d8SPeter Brune 
254fef7b6d8SPeter Brune   Input Parameters:
255fef7b6d8SPeter Brune . snes - the SNES context
256fef7b6d8SPeter Brune 
257fef7b6d8SPeter Brune   Output Parameter:
258fef7b6d8SPeter Brune . outits - number of iterations until termination
259fef7b6d8SPeter Brune 
260fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
261fef7b6d8SPeter Brune */
262b5badacbSBarry Smith static PetscErrorCode SNESSolve_NCG(SNES snes)
263fef7b6d8SPeter Brune {
264dfb256c7SPeter Brune   SNES_NCG             *ncg = (SNES_NCG*)snes->data;
26532cc618eSPeter Brune   Vec                  X,dX,lX,F,dXold;
266bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
26732cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
268fef7b6d8SPeter Brune   PetscInt             maxits, i;
269fef7b6d8SPeter Brune   PetscErrorCode       ierr;
270422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
271f1c6b773SPeter Brune   SNESLineSearch       linesearch;
272b7281c8aSPeter Brune   SNESConvergedReason  reason;
27388d374b2SPeter Brune 
274fef7b6d8SPeter Brune   PetscFunctionBegin;
2756c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
276c579b300SPatrick Farrell 
277fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
278fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
279fef7b6d8SPeter Brune 
280fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
281fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
28232cc618eSPeter Brune   dXold  = snes->work[0];            /* The previous iterate of X */
28388d374b2SPeter Brune   dX     = snes->work[1];            /* the preconditioned direction */
284169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
285fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
286fef7b6d8SPeter Brune 
2877601faf0SJed Brown   ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
2889e764e56SPeter Brune 
289e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
290fef7b6d8SPeter Brune   snes->iter = 0;
291fef7b6d8SPeter Brune   snes->norm = 0.;
292e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
293fef7b6d8SPeter Brune 
294bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
295a71552e2SPeter Brune 
296efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
297be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr);
298efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
299b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
300b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
301b7281c8aSPeter Brune       PetscFunctionReturn(0);
302b7281c8aSPeter Brune     }
303a71552e2SPeter Brune     ierr = VecCopy(dX,F);CHKERRQ(ierr);
304a71552e2SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
305a71552e2SPeter Brune   } else {
306e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
307fef7b6d8SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3081aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
309c1c75074SPeter Brune 
310fef7b6d8SPeter Brune     /* convergence test */
311a71552e2SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
312422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
313a71552e2SPeter Brune     ierr = VecCopy(F,dX);CHKERRQ(ierr);
314a71552e2SPeter Brune   }
315efd4aadfSBarry Smith   if (snes->npc) {
316a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
317be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr);
318efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
319b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
320b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
321b7281c8aSPeter Brune         PetscFunctionReturn(0);
322b7281c8aSPeter Brune       }
323a71552e2SPeter Brune     }
324a71552e2SPeter Brune   }
325a71552e2SPeter Brune   ierr = VecCopy(dX,lX);CHKERRQ(ierr);
32632cc618eSPeter Brune   ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
327a71552e2SPeter Brune 
328a71552e2SPeter Brune 
329e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
330fef7b6d8SPeter Brune   snes->norm = fnorm;
331e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
332a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
333fef7b6d8SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
334fef7b6d8SPeter Brune 
335fef7b6d8SPeter Brune   /* test convergence */
336fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
337fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
338fef7b6d8SPeter Brune 
339fef7b6d8SPeter Brune   /* Call general purpose update function */
340fef7b6d8SPeter Brune   if (snes->ops->update) {
341fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
342fef7b6d8SPeter Brune   }
343fef7b6d8SPeter Brune 
344fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
345bbd5d0b3SPeter Brune 
34609c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
34729c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3480a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
34932cc618eSPeter Brune       ierr = VecCopy(dX, dXold);CHKERRQ(ierr);
35088d374b2SPeter Brune     }
351f1c6b773SPeter Brune     ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr);
352422a814eSBarry Smith     ierr = SNESLineSearchGetReason(linesearch, &lsresult);CHKERRQ(ierr);
353422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
354422a814eSBarry Smith     if (lsresult) {
355fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
356fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3575d10c001SPeter Brune         PetscFunctionReturn(0);
358fef7b6d8SPeter Brune       }
359fef7b6d8SPeter Brune     }
360e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
361fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3625d10c001SPeter Brune       PetscFunctionReturn(0);
363fef7b6d8SPeter Brune     }
364fef7b6d8SPeter Brune     /* Monitor convergence */
365e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
366169e2be8SPeter Brune     snes->iter = i;
367fef7b6d8SPeter Brune     snes->norm = fnorm;
368c1e67a49SFande Kong     snes->xnorm = xnorm;
369c1e67a49SFande Kong     snes->ynorm = ynorm;
370e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
371a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
372fef7b6d8SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
373302dbbaeSPeter Brune 
374fef7b6d8SPeter Brune     /* Test for convergence */
375d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3765d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
377302dbbaeSPeter Brune 
378302dbbaeSPeter Brune     /* Call general purpose update function */
379302dbbaeSPeter Brune     if (snes->ops->update) {
380302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
381302dbbaeSPeter Brune     }
382efd4aadfSBarry Smith     if (snes->npc) {
383a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
384be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr);
385efd4aadfSBarry Smith         ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
386b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
387b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
388b7281c8aSPeter Brune           PetscFunctionReturn(0);
389b7281c8aSPeter Brune         }
390a71552e2SPeter Brune         ierr = VecCopy(dX,F);CHKERRQ(ierr);
391a71552e2SPeter Brune       } else {
392be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr);
393efd4aadfSBarry Smith         ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
394b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
395b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
396b7281c8aSPeter Brune           PetscFunctionReturn(0);
397b7281c8aSPeter Brune         }
398302dbbaeSPeter Brune       }
399a3685310SPeter Brune     } else {
400a3685310SPeter Brune       ierr = VecCopy(F,dX);CHKERRQ(ierr);
401302dbbaeSPeter Brune     }
402302dbbaeSPeter Brune 
40329c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
4040a844d1aSPeter Brune     switch (ncg->type) {
4050a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
40632cc618eSPeter Brune       dXolddotdXold= dXdotdX;
40732cc618eSPeter Brune       ierr         = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
40832cc618eSPeter Brune       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
409dfb256c7SPeter Brune       break;
4100a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
41132cc618eSPeter Brune       dXolddotdXold= dXdotdX;
41232cc618eSPeter Brune       ierr         = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr);
41332cc618eSPeter Brune       ierr         = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr);
41432cc618eSPeter Brune       ierr         = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
41532cc618eSPeter Brune       ierr         = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr);
41632cc618eSPeter Brune       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
417dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
418dfb256c7SPeter Brune       break;
4190a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
42032cc618eSPeter Brune       ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr);
42132cc618eSPeter Brune       ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr);
42232cc618eSPeter Brune       ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr);
42332cc618eSPeter Brune       ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
42432cc618eSPeter Brune       ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
42532cc618eSPeter Brune       ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr);
42632cc618eSPeter Brune       ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr);
42732cc618eSPeter Brune       ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
42832cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
429dfb256c7SPeter Brune       break;
4300a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
43132cc618eSPeter Brune       ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr);
43232cc618eSPeter Brune       ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr);
43332cc618eSPeter Brune       ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
43432cc618eSPeter Brune       ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
43532cc618eSPeter Brune       ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr);
43632cc618eSPeter Brune       ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
43732cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr);
438dfb256c7SPeter Brune       break;
4390a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
44032cc618eSPeter Brune       ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr);
44132cc618eSPeter Brune       ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
44232cc618eSPeter Brune       ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
44332cc618eSPeter Brune       ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
44432cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr);
445dfb256c7SPeter Brune       break;
446dfb256c7SPeter Brune     }
447dfb256c7SPeter Brune     if (ncg->monitor) {
4487904a332SBarry Smith       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);CHKERRQ(ierr);
449dfb256c7SPeter Brune     }
450dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
451fef7b6d8SPeter Brune   }
452fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
453fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
454fef7b6d8SPeter Brune   PetscFunctionReturn(0);
455fef7b6d8SPeter Brune }
456fef7b6d8SPeter Brune 
457fef7b6d8SPeter Brune /*MC
458fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
459fef7b6d8SPeter Brune 
460fef7b6d8SPeter Brune   Level: beginner
461fef7b6d8SPeter Brune 
462fef7b6d8SPeter Brune   Options Database:
4634f02bc6aSBarry Smith +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
46441a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
465b5badacbSBarry Smith -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .
466fef7b6d8SPeter Brune 
46795452b02SPatrick Sanan    Notes:
46895452b02SPatrick Sanan     This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
469fef7b6d8SPeter Brune           gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
470fef7b6d8SPeter Brune           chooses the initial search direction as F(x) for the initial guess x.
471fef7b6d8SPeter Brune 
4722d547940SBarry Smith           Only supports left non-linear preconditioning.
4732d547940SBarry Smith 
474b5badacbSBarry Smith     Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
475b5badacbSBarry Smith     SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR
476b5badacbSBarry Smith 
4774f02bc6aSBarry Smith    References:
47896a0c994SBarry Smith .  1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
4794f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
480cc3c3c0fSMatthew G. Knepley 
481fef7b6d8SPeter Brune 
482b5badacbSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType()
483fef7b6d8SPeter Brune M*/
4848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
485fef7b6d8SPeter Brune {
486fef7b6d8SPeter Brune   PetscErrorCode ierr;
487ea630c6eSPeter Brune   SNES_NCG       * neP;
488fef7b6d8SPeter Brune 
489fef7b6d8SPeter Brune   PetscFunctionBegin;
490fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
491fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
492fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
493fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
494fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
495fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
496fef7b6d8SPeter Brune 
497fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
498efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
499efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
500fef7b6d8SPeter Brune 
5014fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5024fc747eaSLawrence Mitchell 
50388976e71SPeter Brune   if (!snes->tolerancesset) {
5040e444f03SPeter Brune     snes->max_funcs = 30000;
5050e444f03SPeter Brune     snes->max_its   = 10000;
506c522fa08SPeter Brune     snes->stol      = 1e-20;
50788976e71SPeter Brune   }
5080e444f03SPeter Brune 
509b00a9115SJed Brown   ierr         = PetscNewLog(snes,&neP);CHKERRQ(ierr);
510fef7b6d8SPeter Brune   snes->data   = (void*) neP;
5110298fd71SBarry Smith   neP->monitor = NULL;
5120a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
513bdf89e91SBarry Smith   ierr         = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr);
514fef7b6d8SPeter Brune   PetscFunctionReturn(0);
515fef7b6d8SPeter Brune }
516