xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 189a97108d30a3bf7243c9b4a730cb6ba5efe6af)
10a844d1aSPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/
26a6fc655SJed Brown const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0};
3fef7b6d8SPeter Brune 
4fef7b6d8SPeter Brune #undef __FUNCT__
5fef7b6d8SPeter Brune #define __FUNCT__ "SNESReset_NCG"
6fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes)
7fef7b6d8SPeter Brune {
8fef7b6d8SPeter Brune   PetscFunctionBegin;
9fef7b6d8SPeter Brune   PetscFunctionReturn(0);
10fef7b6d8SPeter Brune }
11fef7b6d8SPeter Brune 
121a4f838cSPeter Brune #define SNESLINESEARCHNCGLINEAR "linear"
13bbd5d0b3SPeter Brune 
14fef7b6d8SPeter Brune /*
15fef7b6d8SPeter Brune   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
16fef7b6d8SPeter Brune 
17fef7b6d8SPeter Brune   Input Parameter:
18fef7b6d8SPeter Brune . snes - the SNES context
19fef7b6d8SPeter Brune 
20fef7b6d8SPeter Brune   Application Interface Routine: SNESDestroy()
21fef7b6d8SPeter Brune */
22fef7b6d8SPeter Brune #undef __FUNCT__
23fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG"
24fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes)
25fef7b6d8SPeter Brune {
26fef7b6d8SPeter Brune   PetscErrorCode ierr;
27fef7b6d8SPeter Brune 
28fef7b6d8SPeter Brune   PetscFunctionBegin;
29fef7b6d8SPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
30fef7b6d8SPeter Brune   PetscFunctionReturn(0);
31fef7b6d8SPeter Brune }
32fef7b6d8SPeter Brune 
33fef7b6d8SPeter Brune /*
34fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
35fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
36fef7b6d8SPeter Brune 
37fef7b6d8SPeter Brune    Input Parameters:
38fef7b6d8SPeter Brune +  snes - the SNES context
39fef7b6d8SPeter Brune -  x - the solution vector
40fef7b6d8SPeter Brune 
41fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
42fef7b6d8SPeter Brune  */
43bbd5d0b3SPeter Brune 
445dc0b524SSatish Balay EXTERN_C_BEGIN
45f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch);
465dc0b524SSatish Balay EXTERN_C_END
47bbd5d0b3SPeter Brune 
48fef7b6d8SPeter Brune #undef __FUNCT__
49fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG"
50fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes)
51fef7b6d8SPeter Brune {
52fef7b6d8SPeter Brune   PetscErrorCode ierr;
53fef7b6d8SPeter Brune 
54fef7b6d8SPeter Brune   PetscFunctionBegin;
55bbd5d0b3SPeter Brune   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
56bbd5d0b3SPeter Brune   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
570298fd71SBarry Smith   ierr = SNESLineSearchRegisterDynamic(SNESLINESEARCHNCGLINEAR, NULL,"SNESLineSearchCreate_NCGLinear", SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr);
58fef7b6d8SPeter Brune   PetscFunctionReturn(0);
59fef7b6d8SPeter Brune }
60fef7b6d8SPeter Brune /*
6105b53524SPeter Brune   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
62fef7b6d8SPeter Brune 
63fef7b6d8SPeter Brune   Input Parameter:
64fef7b6d8SPeter Brune . snes - the SNES context
65fef7b6d8SPeter Brune 
66fef7b6d8SPeter Brune   Application Interface Routine: SNESSetFromOptions()
67fef7b6d8SPeter Brune */
68fef7b6d8SPeter Brune #undef __FUNCT__
69fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG"
70fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
71fef7b6d8SPeter Brune {
72dfb256c7SPeter Brune   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
73fef7b6d8SPeter Brune   PetscErrorCode ierr;
740a844d1aSPeter Brune   PetscBool      debug;
75f1c6b773SPeter Brune   SNESLineSearch linesearch;
76a11d90f7SPeter Brune   SNESNCGType    ncgtype=ncg->type;
77fef7b6d8SPeter Brune 
78fef7b6d8SPeter Brune   PetscFunctionBegin;
79fef7b6d8SPeter Brune   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
800298fd71SBarry Smith   ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr);
810298fd71SBarry Smith   ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr);
820a844d1aSPeter Brune   ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr);
83dfb256c7SPeter Brune   if (debug) {
84ce94432eSBarry Smith     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
85dfb256c7SPeter Brune   }
86fef7b6d8SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
879e764e56SPeter Brune   if (!snes->linesearch) {
88f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
899e764e56SPeter Brune     if (!snes->pc) {
901a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
919e764e56SPeter Brune     } else {
921a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
939e764e56SPeter Brune     }
949e764e56SPeter Brune   }
95fef7b6d8SPeter Brune   PetscFunctionReturn(0);
96fef7b6d8SPeter Brune }
97fef7b6d8SPeter Brune 
98fef7b6d8SPeter Brune /*
99fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
100fef7b6d8SPeter Brune 
101fef7b6d8SPeter Brune   Input Parameters:
102fef7b6d8SPeter Brune + SNES - the SNES context
103fef7b6d8SPeter Brune - viewer - visualization context
104fef7b6d8SPeter Brune 
105fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
106fef7b6d8SPeter Brune */
107fef7b6d8SPeter Brune #undef __FUNCT__
108fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
109fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
110fef7b6d8SPeter Brune {
111fef7b6d8SPeter Brune   PetscBool      iascii;
112fef7b6d8SPeter Brune   PetscErrorCode ierr;
113fef7b6d8SPeter Brune 
114fef7b6d8SPeter Brune   PetscFunctionBegin;
115251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
116fef7b6d8SPeter Brune   if (iascii) {
117fef7b6d8SPeter Brune   }
118fef7b6d8SPeter Brune   PetscFunctionReturn(0);
119fef7b6d8SPeter Brune }
120fef7b6d8SPeter Brune 
121fef7b6d8SPeter Brune #undef __FUNCT__
122f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_NCGLinear"
123f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
124ea630c6eSPeter Brune {
125ea630c6eSPeter Brune   PetscScalar    alpha, ptAp;
126bbd5d0b3SPeter Brune   Vec            X, Y, F, W;
127bbd5d0b3SPeter Brune   SNES           snes;
128ea630c6eSPeter Brune   PetscErrorCode ierr;
129bbd5d0b3SPeter Brune   PetscReal      *fnorm, *xnorm, *ynorm;
130ea630c6eSPeter Brune   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
131ea630c6eSPeter Brune 
132ea630c6eSPeter Brune   PetscFunctionBegin;
133f1c6b773SPeter Brune   ierr  = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
134bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
135bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
136bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
137bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
138bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
139bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
140bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
141bbd5d0b3SPeter Brune 
142ea630c6eSPeter Brune   /*
143ea630c6eSPeter Brune 
144d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
145ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
146ea630c6eSPeter Brune    */
147a5892487SPeter Brune   ierr  = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
148ea630c6eSPeter Brune   ierr  = VecDot(F, F, &alpha);CHKERRQ(ierr);
149ea630c6eSPeter Brune   ierr  = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
150ea630c6eSPeter Brune   ierr  = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
151ea630c6eSPeter Brune   alpha = alpha / ptAp;
152ce94432eSBarry Smith   ierr  = PetscPrintf(PetscObjectComm((PetscObject)snes), "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
153bbd5d0b3SPeter Brune   ierr  = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
154bbd5d0b3SPeter Brune   ierr  = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
155bbd5d0b3SPeter Brune 
156bbd5d0b3SPeter Brune   ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
157bbd5d0b3SPeter Brune   ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr);
158bbd5d0b3SPeter Brune   ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr);
159ea630c6eSPeter Brune   PetscFunctionReturn(0);
160ea630c6eSPeter Brune }
161ea630c6eSPeter Brune 
162bbd5d0b3SPeter Brune EXTERN_C_BEGIN
163bbd5d0b3SPeter Brune #undef __FUNCT__
164f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_NCGLinear"
165bbd5d0b3SPeter Brune 
166f1c6b773SPeter Brune PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
167bbd5d0b3SPeter Brune {
168bbd5d0b3SPeter Brune   PetscFunctionBegin;
169f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
1700298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1710298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1720298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1730298fd71SBarry Smith   linesearch->ops->view           = NULL;
1740298fd71SBarry Smith   linesearch->ops->setup          = NULL;
175bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
176bbd5d0b3SPeter Brune }
177bbd5d0b3SPeter Brune EXTERN_C_END
178bbd5d0b3SPeter Brune 
17988d374b2SPeter Brune #undef __FUNCT__
1808c40d5fbSBarry Smith #define __FUNCT__ "SNESNCGComputeYtJtF_Private"
18188d374b2SPeter Brune /*
18288d374b2SPeter Brune 
18388d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
18488d374b2SPeter Brune 
18588d374b2SPeter Brune  */
18667392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
18767392de3SBarry Smith {
18888d374b2SPeter Brune   PetscErrorCode ierr;
18988d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
19067392de3SBarry Smith 
19188d374b2SPeter Brune   PetscFunctionBegin;
19288d374b2SPeter Brune   ierr   = VecDot(F, F, &ftf);CHKERRQ(ierr);
19388d374b2SPeter Brune   ierr   = VecDot(F, Y, &fty);CHKERRQ(ierr);
19488d374b2SPeter Brune   h      = 1e-5*fty / fty;
19588d374b2SPeter Brune   ierr   = VecCopy(X, W);CHKERRQ(ierr);
19688d374b2SPeter Brune   ierr   = VecAXPY(W, -h, Y);CHKERRQ(ierr);          /* this is arbitrary */
19788d374b2SPeter Brune   ierr   = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
19888d374b2SPeter Brune   ierr   = VecDot(G, F, &ftg);CHKERRQ(ierr);
19988d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
20088d374b2SPeter Brune   PetscFunctionReturn(0);
20188d374b2SPeter Brune }
20288d374b2SPeter Brune 
2030a844d1aSPeter Brune #undef __FUNCT__
2040a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType"
2050a844d1aSPeter Brune /*@
2060a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
2070a844d1aSPeter Brune 
2080a844d1aSPeter Brune     Logically Collective on SNES
2090a844d1aSPeter Brune 
2100a844d1aSPeter Brune     Input Parameters:
2110a844d1aSPeter Brune +   snes - the iterative context
2120a844d1aSPeter Brune -   btype - update type
2130a844d1aSPeter Brune 
2140a844d1aSPeter Brune     Options Database:
2150a844d1aSPeter Brune .   -snes_ncg_type<prp,fr,hs,dy,cd>
2160a844d1aSPeter Brune 
2170a844d1aSPeter Brune     Level: intermediate
2180a844d1aSPeter Brune 
2190a844d1aSPeter Brune     SNESNCGSelectTypes:
2200a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2210a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2220a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2230a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2240a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2250a844d1aSPeter Brune 
2260a844d1aSPeter Brune    Notes:
2270a844d1aSPeter Brune    PRP is the default, and the only one that tolerates generalized search directions.
2280a844d1aSPeter Brune 
2290a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set
2300a844d1aSPeter Brune @*/
2310adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
2320adebc6cSBarry Smith {
2330a844d1aSPeter Brune   PetscErrorCode ierr;
2346e111a19SKarl Rupp 
2350a844d1aSPeter Brune   PetscFunctionBegin;
2360a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2370a844d1aSPeter Brune   ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr);
2380a844d1aSPeter Brune   PetscFunctionReturn(0);
2390a844d1aSPeter Brune }
2400a844d1aSPeter Brune 
2410a844d1aSPeter Brune 
2420a844d1aSPeter Brune EXTERN_C_BEGIN
2430a844d1aSPeter Brune #undef __FUNCT__
2440a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType_NCG"
2450adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
2460adebc6cSBarry Smith {
2470a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG*)snes->data;
2486e111a19SKarl Rupp 
2490a844d1aSPeter Brune   PetscFunctionBegin;
2500a844d1aSPeter Brune   ncg->type = btype;
2510a844d1aSPeter Brune   PetscFunctionReturn(0);
2520a844d1aSPeter Brune }
2530a844d1aSPeter Brune EXTERN_C_END
2540a844d1aSPeter Brune 
255fef7b6d8SPeter Brune /*
256fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
257fef7b6d8SPeter Brune 
258fef7b6d8SPeter Brune   Input Parameters:
259fef7b6d8SPeter Brune . snes - the SNES context
260fef7b6d8SPeter Brune 
261fef7b6d8SPeter Brune   Output Parameter:
262fef7b6d8SPeter Brune . outits - number of iterations until termination
263fef7b6d8SPeter Brune 
264fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
265fef7b6d8SPeter Brune */
266fef7b6d8SPeter Brune #undef __FUNCT__
267fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
268fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
269fef7b6d8SPeter Brune {
270dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG*)snes->data;
271bbd5d0b3SPeter Brune   Vec                 X, dX, lX, F, B, Fold;
272bbd5d0b3SPeter Brune   PetscReal           fnorm, ynorm, xnorm, beta = 0.0;
2735d115551SPeter Brune   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
274fef7b6d8SPeter Brune   PetscInt            maxits, i;
275fef7b6d8SPeter Brune   PetscErrorCode      ierr;
276fef7b6d8SPeter Brune   SNESConvergedReason reason;
277fef7b6d8SPeter Brune   PetscBool           lsSuccess = PETSC_TRUE;
278f1c6b773SPeter Brune   SNESLineSearch      linesearch;
27988d374b2SPeter Brune 
280fef7b6d8SPeter Brune   PetscFunctionBegin;
281fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
282fef7b6d8SPeter Brune 
283fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
284fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
2855d115551SPeter Brune   Fold   = snes->work[0];            /* The previous iterate of X */
28688d374b2SPeter Brune   dX     = snes->work[1];            /* the preconditioned direction */
287169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
288fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
289302dbbaeSPeter Brune   B      = snes->vec_rhs;            /* the right hand side */
290fef7b6d8SPeter Brune 
291f1c6b773SPeter Brune   ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
2929e764e56SPeter Brune 
293fef7b6d8SPeter Brune   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
294fef7b6d8SPeter Brune   snes->iter = 0;
295fef7b6d8SPeter Brune   snes->norm = 0.;
296fef7b6d8SPeter Brune   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
297fef7b6d8SPeter Brune 
298bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
299e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
300fef7b6d8SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
301fef7b6d8SPeter Brune     if (snes->domainerror) {
302fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
303fef7b6d8SPeter Brune       PetscFunctionReturn(0);
304fef7b6d8SPeter Brune     }
3051aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
3061aa26658SKarl Rupp 
307e4ed7901SPeter Brune   if (!snes->norm_init_set) {
308fef7b6d8SPeter Brune     /* convergence test */
309fef7b6d8SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
310*189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
311*189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
312*189a9710SBarry Smith       PetscFunctionReturn(0);
313*189a9710SBarry Smith     }
314e4ed7901SPeter Brune   } else {
315e4ed7901SPeter Brune     fnorm               = snes->norm_init;
316e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
317e4ed7901SPeter Brune   }
318fef7b6d8SPeter Brune   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
319fef7b6d8SPeter Brune   snes->norm = fnorm;
320fef7b6d8SPeter Brune   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
321a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
322fef7b6d8SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
323fef7b6d8SPeter Brune 
324fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
325fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
326fef7b6d8SPeter Brune   /* test convergence */
327fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
328fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
329fef7b6d8SPeter Brune 
330fef7b6d8SPeter Brune   /* Call general purpose update function */
331fef7b6d8SPeter Brune   if (snes->ops->update) {
332fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
333fef7b6d8SPeter Brune   }
334fef7b6d8SPeter Brune 
335fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
336bbd5d0b3SPeter Brune 
337c40d0f55SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
33829c94e34SPeter Brune     ierr = VecCopy(X, dX);CHKERRQ(ierr);
339217b9c2eSPeter Brune     ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
340217b9c2eSPeter Brune     ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
34163e7833aSPeter Brune 
34263e7833aSPeter Brune     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,dX,B,0);CHKERRQ(ierr);
34329c94e34SPeter Brune     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
34463e7833aSPeter Brune     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,dX,B,0);CHKERRQ(ierr);
34563e7833aSPeter Brune 
346fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
34751e86f29SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
348fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
349fef7b6d8SPeter Brune       PetscFunctionReturn(0);
350fef7b6d8SPeter Brune     }
35129c94e34SPeter Brune     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
3525d10c001SPeter Brune   } else {
35329c94e34SPeter Brune     ierr = VecCopy(F, dX);CHKERRQ(ierr);
354fef7b6d8SPeter Brune   }
35529c94e34SPeter Brune   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
356bbd5d0b3SPeter Brune   ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
357bbd5d0b3SPeter Brune   /*
35888d374b2SPeter Brune   } else {
3598c40d5fbSBarry Smith     ierr = SNESNCGComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
36088d374b2SPeter Brune   }
361bbd5d0b3SPeter Brune    */
36209c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
363fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
36429c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3650a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
3665d115551SPeter Brune       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
36788d374b2SPeter Brune     }
368f1c6b773SPeter Brune     ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, lX);CHKERRQ(ierr);
369f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr);
370fef7b6d8SPeter Brune     if (!lsSuccess) {
371fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
372fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3735d10c001SPeter Brune         PetscFunctionReturn(0);
374fef7b6d8SPeter Brune       }
375fef7b6d8SPeter Brune     }
376fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
377fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3785d10c001SPeter Brune       PetscFunctionReturn(0);
379fef7b6d8SPeter Brune     }
380fef7b6d8SPeter Brune     if (snes->domainerror) {
381fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
382fef7b6d8SPeter Brune       PetscFunctionReturn(0);
383fef7b6d8SPeter Brune     }
384f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
385fef7b6d8SPeter Brune     /* Monitor convergence */
386fef7b6d8SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
387169e2be8SPeter Brune     snes->iter = i;
388fef7b6d8SPeter Brune     snes->norm = fnorm;
389fef7b6d8SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
390a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
391fef7b6d8SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
392302dbbaeSPeter Brune 
393fef7b6d8SPeter Brune     /* Test for convergence */
394d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3955d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
396302dbbaeSPeter Brune 
397302dbbaeSPeter Brune     /* Call general purpose update function */
398302dbbaeSPeter Brune     if (snes->ops->update) {
399302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
400302dbbaeSPeter Brune     }
401c40d0f55SPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
402302dbbaeSPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
403217b9c2eSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
404217b9c2eSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
405cf949ffaSPeter Brune       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
406302dbbaeSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
40751e86f29SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
408302dbbaeSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
409302dbbaeSPeter Brune         PetscFunctionReturn(0);
410302dbbaeSPeter Brune       }
411eb1825c3SPeter Brune       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
412a3685310SPeter Brune     } else {
413a3685310SPeter Brune       ierr = VecCopy(F, dX);CHKERRQ(ierr);
414302dbbaeSPeter Brune     }
415302dbbaeSPeter Brune 
41629c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
4170a844d1aSPeter Brune     switch (ncg->type) {
4180a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
4195d115551SPeter Brune       dXolddotFold = dXdotF;
420bbd5d0b3SPeter Brune       ierr         = VecDot(dX, dX, &dXdotF);CHKERRQ(ierr);
4215d115551SPeter Brune       beta         = PetscRealPart(dXdotF / dXolddotFold);
422dfb256c7SPeter Brune       break;
4230a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
4245d115551SPeter Brune       dXolddotFold = dXdotF;
42515f0db4aSPeter Brune       ierr         = VecDotBegin(F, dX, &dXdotF);CHKERRQ(ierr);
42615f0db4aSPeter Brune       ierr         = VecDotBegin(Fold, dX, &dXdotFold);CHKERRQ(ierr);
42715f0db4aSPeter Brune       ierr         = VecDotEnd(F, dX, &dXdotF);CHKERRQ(ierr);
42815f0db4aSPeter Brune       ierr         = VecDotEnd(Fold, dX, &dXdotFold);CHKERRQ(ierr);
4295d115551SPeter Brune       beta         = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
430dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
431dfb256c7SPeter Brune       break;
4320a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
43315f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
43415f0db4aSPeter Brune       ierr = VecDotBegin(dX, Fold, &dXdotFold);CHKERRQ(ierr);
43515f0db4aSPeter Brune       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
43615f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
43715f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
43815f0db4aSPeter Brune       ierr = VecDotEnd(dX, Fold, &dXdotFold);CHKERRQ(ierr);
43915f0db4aSPeter Brune       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
44015f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4415d115551SPeter Brune       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
442dfb256c7SPeter Brune       break;
4430a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
44415f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
44515f0db4aSPeter Brune       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
44615f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
44715f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
44815f0db4aSPeter Brune       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
44915f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4505d115551SPeter Brune       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
451dfb256c7SPeter Brune       break;
4520a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
45315f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
45415f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
45515f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
45615f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4575d115551SPeter Brune       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
458dfb256c7SPeter Brune       break;
459dfb256c7SPeter Brune     }
460dfb256c7SPeter Brune     if (ncg->monitor) {
461646217ecSPeter Brune       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
462dfb256c7SPeter Brune     }
463dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
464fef7b6d8SPeter Brune   }
465fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
466fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
467fef7b6d8SPeter Brune   PetscFunctionReturn(0);
468fef7b6d8SPeter Brune }
469fef7b6d8SPeter Brune 
4700a844d1aSPeter Brune 
4710a844d1aSPeter Brune 
472fef7b6d8SPeter Brune /*MC
473fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
474fef7b6d8SPeter Brune 
475fef7b6d8SPeter Brune   Level: beginner
476fef7b6d8SPeter Brune 
477fef7b6d8SPeter Brune   Options Database:
4781867fe5bSPeter Brune +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
47941a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
4801867fe5bSPeter Brune -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
481fef7b6d8SPeter Brune 
482fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
483fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
484fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
485fef7b6d8SPeter Brune 
486fef7b6d8SPeter Brune 
48704d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
488fef7b6d8SPeter Brune M*/
489fef7b6d8SPeter Brune EXTERN_C_BEGIN
490fef7b6d8SPeter Brune #undef __FUNCT__
491fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
492fef7b6d8SPeter Brune PetscErrorCode  SNESCreate_NCG(SNES snes)
493fef7b6d8SPeter Brune {
494fef7b6d8SPeter Brune   PetscErrorCode ierr;
495ea630c6eSPeter Brune   SNES_NCG       * neP;
496fef7b6d8SPeter Brune 
497fef7b6d8SPeter Brune   PetscFunctionBegin;
498fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
499fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
500fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
501fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
502fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
503fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
504fef7b6d8SPeter Brune 
505fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
506fef7b6d8SPeter Brune   snes->usespc  = PETSC_TRUE;
507fef7b6d8SPeter Brune 
50888976e71SPeter Brune   if (!snes->tolerancesset) {
5090e444f03SPeter Brune     snes->max_funcs = 30000;
5100e444f03SPeter Brune     snes->max_its   = 10000;
511c522fa08SPeter Brune     snes->stol      = 1e-20;
51288976e71SPeter Brune   }
5130e444f03SPeter Brune 
514fef7b6d8SPeter Brune   ierr         = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
515fef7b6d8SPeter Brune   snes->data   = (void*) neP;
5160298fd71SBarry Smith   neP->monitor = NULL;
5170a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
5180a844d1aSPeter Brune   ierr         = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNCGSetType_C","SNESNCGSetType_NCG", SNESNCGSetType_NCG);CHKERRQ(ierr);
519fef7b6d8SPeter Brune   PetscFunctionReturn(0);
520fef7b6d8SPeter Brune }
521fef7b6d8SPeter Brune EXTERN_C_END
522