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
SNESDestroy_NCG(SNES snes)4d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NCG(SNES snes)
5d71ae5a4SJacob Faibussowitsch {
6fef7b6d8SPeter Brune PetscFunctionBegin;
72e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL));
89566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data));
93ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10fef7b6d8SPeter Brune }
11fef7b6d8SPeter Brune
SNESSetUp_NCG(SNES snes)12d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NCG(SNES snes)
13d71ae5a4SJacob Faibussowitsch {
14fef7b6d8SPeter Brune PetscFunctionBegin;
159566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 2));
1608401ef6SPierre Jolivet PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
176c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
19fef7b6d8SPeter Brune }
20fef7b6d8SPeter Brune
SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)21d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
22d71ae5a4SJacob Faibussowitsch {
23ea630c6eSPeter Brune PetscScalar alpha, ptAp;
24bbd5d0b3SPeter Brune Vec X, Y, F, W;
25bbd5d0b3SPeter Brune SNES snes;
26bbd5d0b3SPeter Brune PetscReal *fnorm, *xnorm, *ynorm;
27ea630c6eSPeter Brune
28ea630c6eSPeter Brune PetscFunctionBegin;
299566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
30bbd5d0b3SPeter Brune X = linesearch->vec_sol;
31bbd5d0b3SPeter Brune W = linesearch->vec_sol_new;
32bbd5d0b3SPeter Brune F = linesearch->vec_func;
33bbd5d0b3SPeter Brune Y = linesearch->vec_update;
34bbd5d0b3SPeter Brune fnorm = &linesearch->fnorm;
35bbd5d0b3SPeter Brune xnorm = &linesearch->xnorm;
36bbd5d0b3SPeter Brune ynorm = &linesearch->ynorm;
37bbd5d0b3SPeter Brune
3848a46eb9SPierre Jolivet if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));
399105210eSPeter Brune
40ea630c6eSPeter Brune /*
41d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just:
42ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
43ea630c6eSPeter Brune */
449566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
45*76c63389SBarry Smith SNESLineSearchCheckJacobianDomainError(snes, linesearch);
46*76c63389SBarry Smith
479566063dSJacob Faibussowitsch PetscCall(VecDot(F, F, &alpha));
489566063dSJacob Faibussowitsch PetscCall(MatMult(snes->jacobian, Y, W));
499566063dSJacob Faibussowitsch PetscCall(VecDot(Y, W, &ptAp));
50ea630c6eSPeter Brune alpha = alpha / ptAp;
519566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, -alpha, Y));
529566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F));
53bbd5d0b3SPeter Brune
549566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm));
55*76c63389SBarry Smith SNESLineSearchCheckFunctionDomainError(snes, linesearch, *fnorm);
569566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, xnorm));
579566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, ynorm));
583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
59ea630c6eSPeter Brune }
60ea630c6eSPeter Brune
61b5badacbSBarry Smith /*MC
62f6dfbefdSBarry Smith SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG`
63b5badacbSBarry Smith
64b5badacbSBarry 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.
65b5badacbSBarry 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.
66b5badacbSBarry Smith
67420bcc1bSBarry Smith Level: advanced
68420bcc1bSBarry Smith
69f6dfbefdSBarry Smith Notes:
70f6dfbefdSBarry Smith This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian
71b5badacbSBarry Smith
72b5badacbSBarry Smith This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.
73b5badacbSBarry Smith
740241b274SPierre Jolivet .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
75b5badacbSBarry Smith M*/
76b5badacbSBarry Smith
SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)77d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
78d71ae5a4SJacob Faibussowitsch {
79bbd5d0b3SPeter Brune PetscFunctionBegin;
80f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear;
810298fd71SBarry Smith linesearch->ops->destroy = NULL;
820298fd71SBarry Smith linesearch->ops->setfromoptions = NULL;
830298fd71SBarry Smith linesearch->ops->reset = NULL;
840298fd71SBarry Smith linesearch->ops->view = NULL;
850298fd71SBarry Smith linesearch->ops->setup = NULL;
863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
87bbd5d0b3SPeter Brune }
88bbd5d0b3SPeter Brune
SNESSetFromOptions_NCG(SNES snes,PetscOptionItems PetscOptionsObject)89ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems PetscOptionsObject)
90d71ae5a4SJacob Faibussowitsch {
91b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *)snes->data;
92b5badacbSBarry Smith PetscBool debug = PETSC_FALSE;
93b5badacbSBarry Smith SNESNCGType ncgtype = ncg->type;
94b5badacbSBarry Smith SNESLineSearch linesearch;
95b5badacbSBarry Smith
96b5badacbSBarry Smith PetscFunctionBegin;
97d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
989566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
99ad540459SPierre Jolivet if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
1019566063dSJacob Faibussowitsch PetscCall(SNESNCGSetType(snes, ncgtype));
102d0609cedSBarry Smith PetscOptionsHeadEnd();
103b5badacbSBarry Smith if (!snes->linesearch) {
1049566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch));
105ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) {
106b5badacbSBarry Smith if (!snes->npc) {
1079566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
108b5badacbSBarry Smith } else {
109a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT));
110b5badacbSBarry Smith }
111b5badacbSBarry Smith }
112ec786807SJed Brown }
1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
114b5badacbSBarry Smith }
115b5badacbSBarry Smith
SNESView_NCG(SNES snes,PetscViewer viewer)116d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
117d71ae5a4SJacob Faibussowitsch {
118b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *)snes->data;
1199f196a02SMartin Diehl PetscBool isascii;
120b5badacbSBarry Smith
121b5badacbSBarry Smith PetscFunctionBegin;
1229f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1239f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type]));
1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
125b5badacbSBarry Smith }
126b5badacbSBarry Smith
1270a844d1aSPeter Brune /*@
128f6dfbefdSBarry Smith SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`.
1290a844d1aSPeter Brune
130c3339decSBarry Smith Logically Collective
1310a844d1aSPeter Brune
1320a844d1aSPeter Brune Input Parameters:
1330a844d1aSPeter Brune + snes - the iterative context
134ceaaa498SBarry Smith - btype - update type, see `SNESNCGType`
1350a844d1aSPeter Brune
1363c7db156SBarry Smith Options Database Key:
137b5badacbSBarry Smith . -snes_ncg_type <prp,fr,hs,dy,cd> - strategy for selecting algorithm for computing beta
1380a844d1aSPeter Brune
1390a844d1aSPeter Brune Level: intermediate
1400a844d1aSPeter Brune
1410a844d1aSPeter Brune Notes:
142f6dfbefdSBarry Smith `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.
143b5badacbSBarry Smith
144b5badacbSBarry Smith It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
145f6dfbefdSBarry Smith that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`?
1460a844d1aSPeter Brune
147420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
1480a844d1aSPeter Brune @*/
SNESNCGSetType(SNES snes,SNESNCGType btype)149d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
150d71ae5a4SJacob Faibussowitsch {
1510a844d1aSPeter Brune PetscFunctionBegin;
1520a844d1aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
153cac4c232SBarry Smith PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1550a844d1aSPeter Brune }
1560a844d1aSPeter Brune
SNESNCGSetType_NCG(SNES snes,SNESNCGType btype)157d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
158d71ae5a4SJacob Faibussowitsch {
1590a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data;
1606e111a19SKarl Rupp
1610a844d1aSPeter Brune PetscFunctionBegin;
1620a844d1aSPeter Brune ncg->type = btype;
1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1640a844d1aSPeter Brune }
1650a844d1aSPeter Brune
166fef7b6d8SPeter Brune /*
167fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
168fef7b6d8SPeter Brune
1692fe279fdSBarry Smith Input Parameter:
1702fe279fdSBarry Smith . snes - the `SNES` context
171fef7b6d8SPeter Brune
172fef7b6d8SPeter Brune Output Parameter:
173fef7b6d8SPeter Brune . outits - number of iterations until termination
174fef7b6d8SPeter Brune
175fef7b6d8SPeter Brune Application Interface Routine: SNESSolve()
176fef7b6d8SPeter Brune */
SNESSolve_NCG(SNES snes)177d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NCG(SNES snes)
178d71ae5a4SJacob Faibussowitsch {
179dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data;
18032cc618eSPeter Brune Vec X, dX, lX, F, dXold;
181bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0;
18232cc618eSPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
183fef7b6d8SPeter Brune PetscInt maxits, i;
184f1c6b773SPeter Brune SNESLineSearch linesearch;
185b7281c8aSPeter Brune SNESConvergedReason reason;
18688d374b2SPeter Brune
187fef7b6d8SPeter Brune PetscFunctionBegin;
1880b121fc5SBarry 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);
189c579b300SPatrick Farrell
1909566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
191fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING;
192fef7b6d8SPeter Brune
193fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */
194fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */
19532cc618eSPeter Brune dXold = snes->work[0]; /* The previous iterate of X */
19688d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */
197169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */
198fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */
199fef7b6d8SPeter Brune
2009566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch));
2019e764e56SPeter Brune
2029566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
203fef7b6d8SPeter Brune snes->iter = 0;
204fef7b6d8SPeter Brune snes->norm = 0.;
2059566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
206fef7b6d8SPeter Brune
207bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */
208a71552e2SPeter Brune
209efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2109566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2119566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
212b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
213b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
215b7281c8aSPeter Brune }
2169566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, F));
2179566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm));
218a71552e2SPeter Brune } else {
219ac530a7eSPierre Jolivet if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
220ac530a7eSPierre Jolivet else snes->vec_func_init_set = PETSC_FALSE;
221c1c75074SPeter Brune
222fef7b6d8SPeter Brune /* convergence test */
2239566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm));
224*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
2259566063dSJacob Faibussowitsch PetscCall(VecCopy(F, dX));
226a71552e2SPeter Brune }
227efd4aadfSBarry Smith if (snes->npc) {
228a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2299566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, dX));
2309566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
231b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
232b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
234b7281c8aSPeter Brune }
235a71552e2SPeter Brune }
236a71552e2SPeter Brune }
2379566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, lX));
2389566063dSJacob Faibussowitsch PetscCall(VecDot(dX, dX, &dXdotdX));
239a71552e2SPeter Brune
2409566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
241fef7b6d8SPeter Brune snes->norm = fnorm;
2429566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2439566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
244fef7b6d8SPeter Brune
245fef7b6d8SPeter Brune /* test convergence */
2462d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
2472d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm));
2483ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
249fef7b6d8SPeter Brune
250fef7b6d8SPeter Brune /* Call general purpose update function */
251dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter);
252fef7b6d8SPeter Brune
253fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
254bbd5d0b3SPeter Brune
25509c08436SPeter Brune for (i = 1; i < maxits + 1; i++) {
25629c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */
25748a46eb9SPierre Jolivet if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
2589566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
259*76c63389SBarry Smith if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
2609566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
261*76c63389SBarry Smith SNESCheckLineSearchFailure(snes);
262e71169deSBarry Smith if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
263fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
265fef7b6d8SPeter Brune }
266fef7b6d8SPeter Brune /* Monitor convergence */
2679566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
268169e2be8SPeter Brune snes->iter = i;
269fef7b6d8SPeter Brune snes->norm = fnorm;
270c1e67a49SFande Kong snes->xnorm = xnorm;
271c1e67a49SFande Kong snes->ynorm = ynorm;
2729566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2739566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
274302dbbaeSPeter Brune
275fef7b6d8SPeter Brune /* Test for convergence */
2762d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
2772d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2783ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
279302dbbaeSPeter Brune
280302dbbaeSPeter Brune /* Call general purpose update function */
281dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter);
282efd4aadfSBarry Smith if (snes->npc) {
283a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2849566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, dX));
2859566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
286b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
287b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
289b7281c8aSPeter Brune }
2909566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, F));
291a71552e2SPeter Brune } else {
2929566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, dX));
2939566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
294b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
295b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER;
2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
297b7281c8aSPeter Brune }
298302dbbaeSPeter Brune }
299a3685310SPeter Brune } else {
3009566063dSJacob Faibussowitsch PetscCall(VecCopy(F, dX));
301302dbbaeSPeter Brune }
302302dbbaeSPeter Brune
30329c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
3040a844d1aSPeter Brune switch (ncg->type) {
3050a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */
30632cc618eSPeter Brune dXolddotdXold = dXdotdX;
3079566063dSJacob Faibussowitsch PetscCall(VecDot(dX, dX, &dXdotdX));
30832cc618eSPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold);
309dfb256c7SPeter Brune break;
3100a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
31132cc618eSPeter Brune dXolddotdXold = dXdotdX;
312d3981886SPierre Jolivet PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3139566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
3149566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3159566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
316f4f49eeaSPierre Jolivet beta = PetscRealPart((dXdotdX - dXdotdXold) / dXolddotdXold);
317dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */
318dfb256c7SPeter Brune break;
3190a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */
3209566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3219566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
3229566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3239566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3249566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3259566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
3269566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3279566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
32832cc618eSPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
329dfb256c7SPeter Brune break;
3300a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */
3319566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3329566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dX, &lXdotdX));
3339566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3349566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3359566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dX, &lXdotdX));
3369566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
3372f613bf5SBarry Smith beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
338dfb256c7SPeter Brune break;
3390a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */
3409566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX));
3419566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
3429566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX));
3439566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
3442f613bf5SBarry Smith beta = PetscRealPart(dXdotdX / lXdotdXold);
345dfb256c7SPeter Brune break;
346dfb256c7SPeter Brune }
34748a46eb9SPierre Jolivet if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
3489566063dSJacob Faibussowitsch PetscCall(VecAYPX(lX, beta, dX));
349fef7b6d8SPeter Brune }
3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
351fef7b6d8SPeter Brune }
352fef7b6d8SPeter Brune
353fef7b6d8SPeter Brune /*MC
3541d27aa22SBarry Smith SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems {cite}`bruneknepleysmithtu15`.
355fef7b6d8SPeter Brune
356fef7b6d8SPeter Brune Level: beginner
357fef7b6d8SPeter Brune
358f6dfbefdSBarry Smith Options Database Keys:
35962842358SBarry Smith + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is `prp`.
36041a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type.
3611d27aa22SBarry Smith - -snes_ncg_monitor - Print the beta values nonlinear Conjugate-Gradient used in the iteration, .
362fef7b6d8SPeter Brune
36395452b02SPatrick Sanan Notes:
3641d27aa22SBarry Smith This solves the nonlinear system of equations $ F(x) = 0 $ using the nonlinear generalization of the conjugate
365fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
3661d27aa22SBarry Smith chooses the initial search direction as $ F(x) $ for the initial guess $x$.
367fef7b6d8SPeter Brune
3682d547940SBarry Smith Only supports left non-linear preconditioning.
3692d547940SBarry Smith
37062842358SBarry Smith Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with `-npc_snes_type` <type>, `SNESSetNPC()`, or `SNESGetNPC()` then
371a99ef635SJonas Heinzmann `SNESLINESEARCHSECANT` is used. Also supports the special-purpose line search `SNESLINESEARCHNCGLINEAR`
372b5badacbSBarry Smith
373420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
374fef7b6d8SPeter Brune M*/
SNESCreate_NCG(SNES snes)375d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
376d71ae5a4SJacob Faibussowitsch {
377ea630c6eSPeter Brune SNES_NCG *neP;
378fef7b6d8SPeter Brune
379fef7b6d8SPeter Brune PetscFunctionBegin;
380fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG;
381fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG;
382fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG;
383fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG;
384fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG;
385fef7b6d8SPeter Brune
386fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE;
387efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE;
388efd4aadfSBarry Smith snes->npcside = PC_LEFT;
389fef7b6d8SPeter Brune
3904fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE;
3914fc747eaSLawrence Mitchell
39277e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes));
39377e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_funcs, 30000);
39477e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_its, 10000);
39577e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, stol, 1e-20);
3960e444f03SPeter Brune
3974dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP));
398fef7b6d8SPeter Brune snes->data = (void *)neP;
3990298fd71SBarry Smith neP->monitor = NULL;
4000a844d1aSPeter Brune neP->type = SNES_NCG_PRP;
4019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
403fef7b6d8SPeter Brune }
404