xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 0adebc6ce4218bb5a986f93c65179d861410e0be)
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 
9fef7b6d8SPeter Brune   PetscFunctionBegin;
10fef7b6d8SPeter Brune   PetscFunctionReturn(0);
11fef7b6d8SPeter Brune }
12fef7b6d8SPeter Brune 
131a4f838cSPeter Brune #define SNESLINESEARCHNCGLINEAR "linear"
14bbd5d0b3SPeter Brune 
15fef7b6d8SPeter Brune /*
16fef7b6d8SPeter Brune   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
17fef7b6d8SPeter Brune 
18fef7b6d8SPeter Brune   Input Parameter:
19fef7b6d8SPeter Brune . snes - the SNES context
20fef7b6d8SPeter Brune 
21fef7b6d8SPeter Brune   Application Interface Routine: SNESDestroy()
22fef7b6d8SPeter Brune */
23fef7b6d8SPeter Brune #undef __FUNCT__
24fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG"
25fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes)
26fef7b6d8SPeter Brune {
27fef7b6d8SPeter Brune   PetscErrorCode   ierr;
28fef7b6d8SPeter Brune 
29fef7b6d8SPeter Brune   PetscFunctionBegin;
30fef7b6d8SPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
31fef7b6d8SPeter Brune   PetscFunctionReturn(0);
32fef7b6d8SPeter Brune }
33fef7b6d8SPeter Brune 
34fef7b6d8SPeter Brune /*
35fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
36fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
37fef7b6d8SPeter Brune 
38fef7b6d8SPeter Brune    Input Parameters:
39fef7b6d8SPeter Brune +  snes - the SNES context
40fef7b6d8SPeter Brune -  x - the solution vector
41fef7b6d8SPeter Brune 
42fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
43fef7b6d8SPeter Brune  */
44bbd5d0b3SPeter Brune 
455dc0b524SSatish Balay EXTERN_C_BEGIN
46f1c6b773SPeter Brune extern PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch);
475dc0b524SSatish Balay EXTERN_C_END
48bbd5d0b3SPeter Brune 
49fef7b6d8SPeter Brune #undef __FUNCT__
50fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG"
51fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes)
52fef7b6d8SPeter Brune {
53fef7b6d8SPeter Brune   PetscErrorCode ierr;
54fef7b6d8SPeter Brune 
55fef7b6d8SPeter Brune   PetscFunctionBegin;
56bbd5d0b3SPeter Brune   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
57bbd5d0b3SPeter Brune   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
581a4f838cSPeter Brune   ierr = SNESLineSearchRegisterDynamic(SNESLINESEARCHNCGLINEAR, PETSC_NULL,"SNESLineSearchCreate_NCGLinear", SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr);
59bbd5d0b3SPeter Brune 
60fef7b6d8SPeter Brune   PetscFunctionReturn(0);
61fef7b6d8SPeter Brune }
62fef7b6d8SPeter Brune /*
6305b53524SPeter Brune   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
64fef7b6d8SPeter Brune 
65fef7b6d8SPeter Brune   Input Parameter:
66fef7b6d8SPeter Brune . snes - the SNES context
67fef7b6d8SPeter Brune 
68fef7b6d8SPeter Brune   Application Interface Routine: SNESSetFromOptions()
69fef7b6d8SPeter Brune */
70fef7b6d8SPeter Brune #undef __FUNCT__
71fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG"
72fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
73fef7b6d8SPeter Brune {
74dfb256c7SPeter Brune   SNES_NCG           *ncg     = (SNES_NCG *)snes->data;
75fef7b6d8SPeter Brune   PetscErrorCode     ierr;
760a844d1aSPeter Brune   PetscBool          debug;
77f1c6b773SPeter Brune   SNESLineSearch     linesearch;
78a11d90f7SPeter Brune   SNESNCGType        ncgtype=ncg->type;
79fef7b6d8SPeter Brune 
80fef7b6d8SPeter Brune   PetscFunctionBegin;
81fef7b6d8SPeter Brune   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
82dfb256c7SPeter Brune   ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
830a844d1aSPeter Brune   ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,PETSC_NULL);CHKERRQ(ierr);
840a844d1aSPeter Brune   ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr);
85dfb256c7SPeter Brune   if (debug) {
86dfb256c7SPeter Brune     ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
87dfb256c7SPeter Brune   }
88fef7b6d8SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
899e764e56SPeter Brune   if (!snes->linesearch) {
90f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
919e764e56SPeter Brune     if (!snes->pc) {
921a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr);
939e764e56SPeter Brune     } else {
941a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
959e764e56SPeter Brune     }
969e764e56SPeter Brune   }
97fef7b6d8SPeter Brune   PetscFunctionReturn(0);
98fef7b6d8SPeter Brune }
99fef7b6d8SPeter Brune 
100fef7b6d8SPeter Brune /*
101fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
102fef7b6d8SPeter Brune 
103fef7b6d8SPeter Brune   Input Parameters:
104fef7b6d8SPeter Brune + SNES - the SNES context
105fef7b6d8SPeter Brune - viewer - visualization context
106fef7b6d8SPeter Brune 
107fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
108fef7b6d8SPeter Brune */
109fef7b6d8SPeter Brune #undef __FUNCT__
110fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
111fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
112fef7b6d8SPeter Brune {
113fef7b6d8SPeter Brune   PetscBool        iascii;
114fef7b6d8SPeter Brune   PetscErrorCode   ierr;
115fef7b6d8SPeter Brune 
116fef7b6d8SPeter Brune   PetscFunctionBegin;
117251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
118fef7b6d8SPeter Brune   if (iascii) {
119fef7b6d8SPeter Brune   }
120fef7b6d8SPeter Brune   PetscFunctionReturn(0);
121fef7b6d8SPeter Brune }
122fef7b6d8SPeter Brune 
123fef7b6d8SPeter Brune #undef __FUNCT__
124f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_NCGLinear"
125f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
126ea630c6eSPeter Brune {
127ea630c6eSPeter Brune   PetscScalar      alpha, ptAp;
128bbd5d0b3SPeter Brune   Vec              X, Y, F, W;
129bbd5d0b3SPeter Brune   SNES             snes;
130ea630c6eSPeter Brune   PetscErrorCode   ierr;
131bbd5d0b3SPeter Brune   PetscReal        *fnorm, *xnorm, *ynorm;
132ea630c6eSPeter Brune   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
133ea630c6eSPeter Brune 
134ea630c6eSPeter Brune   PetscFunctionBegin;
135bbd5d0b3SPeter Brune 
136f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
137bbd5d0b3SPeter Brune   X = linesearch->vec_sol;
138bbd5d0b3SPeter Brune   W = linesearch->vec_sol_new;
139bbd5d0b3SPeter Brune   F = linesearch->vec_func;
140bbd5d0b3SPeter Brune   Y = linesearch->vec_update;
141bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
142bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
143bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
144bbd5d0b3SPeter Brune 
145ea630c6eSPeter Brune   /*
146ea630c6eSPeter Brune 
147d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
148ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
149ea630c6eSPeter Brune    */
150a5892487SPeter Brune   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
151ea630c6eSPeter Brune   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
152ea630c6eSPeter Brune   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
153ea630c6eSPeter Brune   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
154ea630c6eSPeter Brune   alpha = alpha / ptAp;
155d9fe6452SJed Brown   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
156bbd5d0b3SPeter Brune   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
157bbd5d0b3SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
158bbd5d0b3SPeter Brune 
159bbd5d0b3SPeter Brune   ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
160bbd5d0b3SPeter Brune   ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr);
161bbd5d0b3SPeter Brune   ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr);
162ea630c6eSPeter Brune  PetscFunctionReturn(0);
163ea630c6eSPeter Brune }
164ea630c6eSPeter Brune 
165bbd5d0b3SPeter Brune EXTERN_C_BEGIN
166bbd5d0b3SPeter Brune #undef __FUNCT__
167f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_NCGLinear"
168bbd5d0b3SPeter Brune 
169f1c6b773SPeter Brune PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
170bbd5d0b3SPeter Brune {
171bbd5d0b3SPeter Brune   PetscFunctionBegin;
172f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
173bbd5d0b3SPeter Brune   linesearch->ops->destroy        = PETSC_NULL;
174bbd5d0b3SPeter Brune   linesearch->ops->setfromoptions = PETSC_NULL;
175bbd5d0b3SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
176bbd5d0b3SPeter Brune   linesearch->ops->view           = PETSC_NULL;
177bbd5d0b3SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
178bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
179bbd5d0b3SPeter Brune }
180bbd5d0b3SPeter Brune EXTERN_C_END
181bbd5d0b3SPeter Brune 
18288d374b2SPeter Brune #undef __FUNCT__
1838c40d5fbSBarry Smith #define __FUNCT__ "SNESNCGComputeYtJtF_Private"
18488d374b2SPeter Brune /*
18588d374b2SPeter Brune 
18688d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
18788d374b2SPeter Brune 
18888d374b2SPeter Brune  */
18967392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
19067392de3SBarry Smith  {
19188d374b2SPeter Brune   PetscErrorCode ierr;
19288d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
19367392de3SBarry Smith 
19488d374b2SPeter Brune   PetscFunctionBegin;
19588d374b2SPeter Brune   ierr = VecDot(F, F, &ftf);CHKERRQ(ierr);
19688d374b2SPeter Brune   ierr = VecDot(F, Y, &fty);CHKERRQ(ierr);
19788d374b2SPeter Brune   h = 1e-5*fty / fty;
19888d374b2SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
19988d374b2SPeter Brune   ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr);            /* this is arbitrary */
20088d374b2SPeter Brune   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
20188d374b2SPeter Brune   ierr = VecDot(G, F, &ftg);CHKERRQ(ierr);
20288d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
20388d374b2SPeter Brune   PetscFunctionReturn(0);
20488d374b2SPeter Brune }
20588d374b2SPeter Brune 
2060a844d1aSPeter Brune #undef __FUNCT__
2070a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType"
2080a844d1aSPeter Brune /*@
2090a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
2100a844d1aSPeter Brune 
2110a844d1aSPeter Brune     Logically Collective on SNES
2120a844d1aSPeter Brune 
2130a844d1aSPeter Brune     Input Parameters:
2140a844d1aSPeter Brune +   snes - the iterative context
2150a844d1aSPeter Brune -   btype - update type
2160a844d1aSPeter Brune 
2170a844d1aSPeter Brune     Options Database:
2180a844d1aSPeter Brune .   -snes_ncg_type<prp,fr,hs,dy,cd>
2190a844d1aSPeter Brune 
2200a844d1aSPeter Brune     Level: intermediate
2210a844d1aSPeter Brune 
2220a844d1aSPeter Brune     SNESNCGSelectTypes:
2230a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2240a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2250a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2260a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2270a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2280a844d1aSPeter Brune 
2290a844d1aSPeter Brune    Notes:
2300a844d1aSPeter Brune    PRP is the default, and the only one that tolerates generalized search directions.
2310a844d1aSPeter Brune 
2320a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set
2330a844d1aSPeter Brune @*/
234*0adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
235*0adebc6cSBarry Smith {
2360a844d1aSPeter Brune   PetscErrorCode ierr;
2370a844d1aSPeter Brune   PetscFunctionBegin;
2380a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2390a844d1aSPeter Brune   ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr);
2400a844d1aSPeter Brune   PetscFunctionReturn(0);
2410a844d1aSPeter Brune }
2420a844d1aSPeter Brune 
2430a844d1aSPeter Brune 
2440a844d1aSPeter Brune EXTERN_C_BEGIN
2450a844d1aSPeter Brune #undef __FUNCT__
2460a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType_NCG"
247*0adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
248*0adebc6cSBarry Smith {
2490a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG *)snes->data;
2500a844d1aSPeter Brune   PetscFunctionBegin;
2510a844d1aSPeter Brune   ncg->type = btype;
2520a844d1aSPeter Brune   PetscFunctionReturn(0);
2530a844d1aSPeter Brune }
2540a844d1aSPeter Brune EXTERN_C_END
2550a844d1aSPeter Brune 
256fef7b6d8SPeter Brune /*
257fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
258fef7b6d8SPeter Brune 
259fef7b6d8SPeter Brune   Input Parameters:
260fef7b6d8SPeter Brune . snes - the SNES context
261fef7b6d8SPeter Brune 
262fef7b6d8SPeter Brune   Output Parameter:
263fef7b6d8SPeter Brune . outits - number of iterations until termination
264fef7b6d8SPeter Brune 
265fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
266fef7b6d8SPeter Brune */
267fef7b6d8SPeter Brune #undef __FUNCT__
268fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
269fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
270fef7b6d8SPeter Brune {
271dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
272bbd5d0b3SPeter Brune   Vec                 X, dX, lX, F, B, Fold;
273bbd5d0b3SPeter Brune   PetscReal           fnorm, ynorm, xnorm, beta = 0.0;
2745d115551SPeter Brune   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
275fef7b6d8SPeter Brune   PetscInt            maxits, i;
276fef7b6d8SPeter Brune   PetscErrorCode      ierr;
277fef7b6d8SPeter Brune   SNESConvergedReason reason;
278fef7b6d8SPeter Brune   PetscBool           lsSuccess = PETSC_TRUE;
279f1c6b773SPeter Brune   SNESLineSearch     linesearch;
28088d374b2SPeter Brune 
281fef7b6d8SPeter Brune   PetscFunctionBegin;
282fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
283fef7b6d8SPeter Brune 
284fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
285fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
2865d115551SPeter Brune   Fold   = snes->work[0];            /* The previous iterate of X */
28788d374b2SPeter Brune   dX     = snes->work[1];            /* the preconditioned direction */
288169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
289fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
290302dbbaeSPeter Brune   B      = snes->vec_rhs;            /* the right hand side */
291fef7b6d8SPeter Brune 
292f1c6b773SPeter Brune   ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
2939e764e56SPeter Brune 
294fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
295fef7b6d8SPeter Brune   snes->iter = 0;
296fef7b6d8SPeter Brune   snes->norm = 0.;
297fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
298fef7b6d8SPeter Brune 
299bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
300e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
301fef7b6d8SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
302fef7b6d8SPeter Brune     if (snes->domainerror) {
303fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
304fef7b6d8SPeter Brune       PetscFunctionReturn(0);
305fef7b6d8SPeter Brune     }
306e4ed7901SPeter Brune   } else {
307e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
308e4ed7901SPeter Brune   }
309e4ed7901SPeter Brune   if (!snes->norm_init_set) {
310fef7b6d8SPeter Brune   /* convergence test */
311fef7b6d8SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
312fef7b6d8SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
313e4ed7901SPeter Brune   } else {
314e4ed7901SPeter Brune     fnorm = snes->norm_init;
315e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
316e4ed7901SPeter Brune   }
317fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
318fef7b6d8SPeter Brune   snes->norm = fnorm;
319fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
320fef7b6d8SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
321fef7b6d8SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
322fef7b6d8SPeter Brune 
323fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
324fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
325fef7b6d8SPeter Brune   /* test convergence */
326fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
327fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
328fef7b6d8SPeter Brune 
329fef7b6d8SPeter Brune   /* Call general purpose update function */
330fef7b6d8SPeter Brune   if (snes->ops->update) {
331fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
332fef7b6d8SPeter Brune  }
333fef7b6d8SPeter Brune 
334fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
335bbd5d0b3SPeter Brune 
336c40d0f55SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
33729c94e34SPeter Brune     ierr = VecCopy(X, dX);CHKERRQ(ierr);
338217b9c2eSPeter Brune     ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
339217b9c2eSPeter Brune     ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
34029c94e34SPeter Brune     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
341fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
34251e86f29SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
343fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
344fef7b6d8SPeter Brune       PetscFunctionReturn(0);
345fef7b6d8SPeter Brune     }
34629c94e34SPeter Brune     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
3475d10c001SPeter Brune   } else {
34829c94e34SPeter Brune     ierr = VecCopy(F, dX);CHKERRQ(ierr);
349fef7b6d8SPeter Brune   }
35029c94e34SPeter Brune   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
351bbd5d0b3SPeter Brune   ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
352bbd5d0b3SPeter Brune   /*
35388d374b2SPeter Brune   } else {
3548c40d5fbSBarry Smith     ierr = SNESNCGComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
35588d374b2SPeter Brune   }
356bbd5d0b3SPeter Brune    */
35709c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
358fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
35929c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3600a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
3615d115551SPeter Brune       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
36288d374b2SPeter Brune     }
363f1c6b773SPeter Brune     ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, lX);CHKERRQ(ierr);
364f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr);
365fef7b6d8SPeter Brune     if (!lsSuccess) {
366fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
367fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3685d10c001SPeter Brune         PetscFunctionReturn(0);
369fef7b6d8SPeter Brune       }
370fef7b6d8SPeter Brune     }
371fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
372fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3735d10c001SPeter Brune       PetscFunctionReturn(0);
374fef7b6d8SPeter Brune     }
375fef7b6d8SPeter Brune     if (snes->domainerror) {
376fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
377fef7b6d8SPeter Brune       PetscFunctionReturn(0);
378fef7b6d8SPeter Brune     }
379f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
380fef7b6d8SPeter Brune     /* Monitor convergence */
381fef7b6d8SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
382169e2be8SPeter Brune     snes->iter = i;
383fef7b6d8SPeter Brune     snes->norm = fnorm;
384fef7b6d8SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
385fef7b6d8SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
386fef7b6d8SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
387302dbbaeSPeter Brune 
388fef7b6d8SPeter Brune     /* Test for convergence */
389d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3905d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
391302dbbaeSPeter Brune 
392302dbbaeSPeter Brune     /* Call general purpose update function */
393302dbbaeSPeter Brune     if (snes->ops->update) {
394302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
395302dbbaeSPeter Brune     }
396c40d0f55SPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
397302dbbaeSPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
398217b9c2eSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
399217b9c2eSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
400cf949ffaSPeter Brune       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
401302dbbaeSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
40251e86f29SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
403302dbbaeSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
404302dbbaeSPeter Brune         PetscFunctionReturn(0);
405302dbbaeSPeter Brune       }
406eb1825c3SPeter Brune       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
407a3685310SPeter Brune     } else {
408a3685310SPeter Brune       ierr = VecCopy(F, dX);CHKERRQ(ierr);
409302dbbaeSPeter Brune     }
410302dbbaeSPeter Brune 
41129c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
4120a844d1aSPeter Brune     switch(ncg->type) {
4130a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
4145d115551SPeter Brune       dXolddotFold = dXdotF;
415bbd5d0b3SPeter Brune       ierr = VecDot(dX, dX, &dXdotF);CHKERRQ(ierr);
4165d115551SPeter Brune       beta = PetscRealPart(dXdotF / dXolddotFold);
417dfb256c7SPeter Brune       break;
4180a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
4195d115551SPeter Brune       dXolddotFold = dXdotF;
42015f0db4aSPeter Brune       ierr = VecDotBegin(F, dX, &dXdotF);CHKERRQ(ierr);
42115f0db4aSPeter Brune       ierr = VecDotBegin(Fold, dX, &dXdotFold);CHKERRQ(ierr);
42215f0db4aSPeter Brune       ierr = VecDotEnd(F, dX, &dXdotF);CHKERRQ(ierr);
42315f0db4aSPeter Brune       ierr = VecDotEnd(Fold, dX, &dXdotFold);CHKERRQ(ierr);
4245d115551SPeter Brune       beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
425dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
426dfb256c7SPeter Brune       break;
4270a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
42815f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
42915f0db4aSPeter Brune       ierr = VecDotBegin(dX, Fold, &dXdotFold);CHKERRQ(ierr);
43015f0db4aSPeter Brune       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
43115f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
43215f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
43315f0db4aSPeter Brune       ierr = VecDotEnd(dX, Fold, &dXdotFold);CHKERRQ(ierr);
43415f0db4aSPeter Brune       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
43515f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4365d115551SPeter Brune       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
437dfb256c7SPeter Brune       break;
4380a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
43915f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
44015f0db4aSPeter Brune       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
44115f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
44215f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
44315f0db4aSPeter Brune       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
44415f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4455d115551SPeter Brune       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
446dfb256c7SPeter Brune       break;
4470a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
44815f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
44915f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
45015f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
45115f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4525d115551SPeter Brune       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
453dfb256c7SPeter Brune       break;
454dfb256c7SPeter Brune     }
455dfb256c7SPeter Brune     if (ncg->monitor) {
456646217ecSPeter Brune       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
457dfb256c7SPeter Brune     }
458dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
459fef7b6d8SPeter Brune   }
460fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
461fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
462fef7b6d8SPeter Brune   PetscFunctionReturn(0);
463fef7b6d8SPeter Brune }
464fef7b6d8SPeter Brune 
4650a844d1aSPeter Brune 
4660a844d1aSPeter Brune 
467fef7b6d8SPeter Brune /*MC
468fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
469fef7b6d8SPeter Brune 
470fef7b6d8SPeter Brune   Level: beginner
471fef7b6d8SPeter Brune 
472fef7b6d8SPeter Brune   Options Database:
4731867fe5bSPeter Brune +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
47441a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
4751867fe5bSPeter Brune -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
476fef7b6d8SPeter Brune 
477fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
478fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
479fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
480fef7b6d8SPeter Brune 
481fef7b6d8SPeter Brune 
48204d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
483fef7b6d8SPeter Brune M*/
484fef7b6d8SPeter Brune EXTERN_C_BEGIN
485fef7b6d8SPeter Brune #undef __FUNCT__
486fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
487fef7b6d8SPeter Brune PetscErrorCode  SNESCreate_NCG(SNES snes)
488fef7b6d8SPeter Brune {
489fef7b6d8SPeter Brune   PetscErrorCode   ierr;
490ea630c6eSPeter Brune   SNES_NCG * neP;
491fef7b6d8SPeter Brune 
492fef7b6d8SPeter Brune   PetscFunctionBegin;
493fef7b6d8SPeter Brune   snes->ops->destroy         = SNESDestroy_NCG;
494fef7b6d8SPeter Brune   snes->ops->setup           = SNESSetUp_NCG;
495fef7b6d8SPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
496fef7b6d8SPeter Brune   snes->ops->view            = SNESView_NCG;
497fef7b6d8SPeter Brune   snes->ops->solve           = SNESSolve_NCG;
498fef7b6d8SPeter Brune   snes->ops->reset           = SNESReset_NCG;
499fef7b6d8SPeter Brune 
500fef7b6d8SPeter Brune   snes->usesksp              = PETSC_FALSE;
501fef7b6d8SPeter Brune   snes->usespc               = PETSC_TRUE;
502fef7b6d8SPeter Brune 
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 
509fef7b6d8SPeter Brune   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
510fef7b6d8SPeter Brune   snes->data = (void*) neP;
511dfb256c7SPeter Brune   neP->monitor = PETSC_NULL;
5120a844d1aSPeter Brune   neP->type = SNES_NCG_PRP;
5130a844d1aSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNCGSetType_C","SNESNCGSetType_NCG", SNESNCGSetType_NCG);CHKERRQ(ierr);
514fef7b6d8SPeter Brune   PetscFunctionReturn(0);
515fef7b6d8SPeter Brune }
516fef7b6d8SPeter Brune EXTERN_C_END
517