xref: /petsc/src/snes/impls/ncg/snesncg.c (revision cf9edfecb14be6ea29b532008cdf7dfb760d0bf7)
1fef7b6d8SPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h>
2fef7b6d8SPeter Brune 
3fef7b6d8SPeter Brune #undef __FUNCT__
4fef7b6d8SPeter Brune #define __FUNCT__ "SNESReset_NCG"
5fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes)
6fef7b6d8SPeter Brune {
7fef7b6d8SPeter Brune   PetscErrorCode ierr;
8fef7b6d8SPeter Brune 
9fef7b6d8SPeter Brune   PetscFunctionBegin;
10fef7b6d8SPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
11fef7b6d8SPeter Brune   PetscFunctionReturn(0);
12fef7b6d8SPeter Brune }
13fef7b6d8SPeter 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 = SNESReset_NCG(snes);CHKERRQ(ierr);
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  */
44fef7b6d8SPeter Brune #undef __FUNCT__
45fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG"
46fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes)
47fef7b6d8SPeter Brune {
48fef7b6d8SPeter Brune   PetscErrorCode ierr;
49fef7b6d8SPeter Brune 
50fef7b6d8SPeter Brune   PetscFunctionBegin;
51a5892487SPeter Brune   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
52fef7b6d8SPeter Brune   PetscFunctionReturn(0);
53fef7b6d8SPeter Brune }
54fef7b6d8SPeter Brune /*
55fef7b6d8SPeter Brune   SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method.
56fef7b6d8SPeter Brune 
57fef7b6d8SPeter Brune   Input Parameter:
58fef7b6d8SPeter Brune . snes - the SNES context
59fef7b6d8SPeter Brune 
60fef7b6d8SPeter Brune   Application Interface Routine: SNESSetFromOptions()
61fef7b6d8SPeter Brune */
62fef7b6d8SPeter Brune #undef __FUNCT__
63fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG"
64fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
65fef7b6d8SPeter Brune {
66fef7b6d8SPeter Brune   PetscErrorCode     ierr;
67fef7b6d8SPeter Brune 
68fef7b6d8SPeter Brune   PetscFunctionBegin;
69fef7b6d8SPeter Brune   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
70fef7b6d8SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
71fef7b6d8SPeter Brune   PetscFunctionReturn(0);
72fef7b6d8SPeter Brune }
73fef7b6d8SPeter Brune 
74fef7b6d8SPeter Brune /*
75fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
76fef7b6d8SPeter Brune 
77fef7b6d8SPeter Brune   Input Parameters:
78fef7b6d8SPeter Brune + SNES - the SNES context
79fef7b6d8SPeter Brune - viewer - visualization context
80fef7b6d8SPeter Brune 
81fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
82fef7b6d8SPeter Brune */
83fef7b6d8SPeter Brune #undef __FUNCT__
84fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
85fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
86fef7b6d8SPeter Brune {
87fef7b6d8SPeter Brune   PetscBool        iascii;
88fef7b6d8SPeter Brune   PetscErrorCode   ierr;
89fef7b6d8SPeter Brune 
90fef7b6d8SPeter Brune   PetscFunctionBegin;
91fef7b6d8SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
92fef7b6d8SPeter Brune   if (iascii) {
93ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
94fef7b6d8SPeter Brune   }
95fef7b6d8SPeter Brune   PetscFunctionReturn(0);
96fef7b6d8SPeter Brune }
97fef7b6d8SPeter Brune 
98fef7b6d8SPeter Brune #undef __FUNCT__
99fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_NCG"
100fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchQuadratic_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
101fef7b6d8SPeter Brune {
102fef7b6d8SPeter Brune   PetscInt         i;
103fef7b6d8SPeter Brune   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
104fef7b6d8SPeter Brune   PetscReal        norms[3];
105fef7b6d8SPeter Brune   PetscReal        alpha,a,b;
106fef7b6d8SPeter Brune   PetscErrorCode   ierr;
107fef7b6d8SPeter Brune 
108fef7b6d8SPeter Brune   PetscFunctionBegin;
109fef7b6d8SPeter Brune   norms[0]  = fnorm;
110fef7b6d8SPeter Brune   for(i=1; i < 3; ++i) {
111fef7b6d8SPeter Brune     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
112a5892487SPeter Brune     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
113a5892487SPeter Brune     ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr);
114fef7b6d8SPeter Brune   }
115fef7b6d8SPeter Brune   for(i = 0; i < 3; ++i) {
116fef7b6d8SPeter Brune     norms[i] = PetscSqr(norms[i]);
117fef7b6d8SPeter Brune   }
118fef7b6d8SPeter Brune   /* Fit a quadratic:
119fef7b6d8SPeter Brune        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
120fef7b6d8SPeter Brune        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
121fef7b6d8SPeter Brune        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
122fef7b6d8SPeter Brune        c = y_0
123fef7b6d8SPeter Brune        x_min = -b/2a
124fef7b6d8SPeter Brune 
125fef7b6d8SPeter Brune        If we let x_{0,1,2} = 0, 0.5, 1.0
126fef7b6d8SPeter Brune        a = 2 y_2 - 4 y_1 + 2 y_0
127fef7b6d8SPeter Brune        b =  -y_2 + 4 y_1 - 3 y_0
128fef7b6d8SPeter Brune        c =   y_0
129fef7b6d8SPeter Brune   */
130fef7b6d8SPeter Brune   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
131fef7b6d8SPeter Brune   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
132fef7b6d8SPeter Brune   /* Check for positive a (concave up) */
133fef7b6d8SPeter Brune   if (a >= 0.0) {
134fef7b6d8SPeter Brune     alpha = -b/(2.0*a);
135fef7b6d8SPeter Brune     alpha = PetscMin(alpha, alphas[2]);
136fef7b6d8SPeter Brune     alpha = PetscMax(alpha, alphas[0]);
137fef7b6d8SPeter Brune   } else {
138fef7b6d8SPeter Brune     alpha = 1.0;
139fef7b6d8SPeter Brune   }
140ea630c6eSPeter Brune   if (snes->ls_monitor) {
141ea630c6eSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
142ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
143ea630c6eSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
144fef7b6d8SPeter Brune   }
145a5892487SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
146a5892487SPeter Brune   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
147fef7b6d8SPeter Brune   if (alpha != 1.0) {
148a5892487SPeter Brune     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
149a5892487SPeter Brune     ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
150fef7b6d8SPeter Brune   } else {
151fef7b6d8SPeter Brune     *gnorm = PetscSqrtReal(norms[2]);
152fef7b6d8SPeter Brune   }
153fef7b6d8SPeter Brune   if (alpha == 0.0) *flag = PETSC_FALSE;
154fef7b6d8SPeter Brune   else              *flag = PETSC_TRUE;
155fef7b6d8SPeter Brune   PetscFunctionReturn(0);
156fef7b6d8SPeter Brune }
157fef7b6d8SPeter Brune 
158ea630c6eSPeter Brune 
159ea630c6eSPeter Brune 
160ea630c6eSPeter Brune #undef __FUNCT__
161ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchExactLinear_NCG"
162ea630c6eSPeter Brune PetscErrorCode SNESLineSearchExactLinear_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
163ea630c6eSPeter Brune {
164ea630c6eSPeter Brune   PetscScalar      alpha, ptAp;
165ea630c6eSPeter Brune   PetscErrorCode   ierr;
166ea630c6eSPeter Brune   /* SNES_NCG *ncg =  (SNES_NCG *) snes->data; */
167ea630c6eSPeter Brune   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
168ea630c6eSPeter Brune 
169ea630c6eSPeter Brune   PetscFunctionBegin;
170ea630c6eSPeter Brune   /*
171ea630c6eSPeter Brune 
172ea630c6eSPeter Brune    The exact step size for linear CG is just:
173ea630c6eSPeter Brune 
174ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
175ea630c6eSPeter Brune 
176ea630c6eSPeter Brune    */
177ea630c6eSPeter Brune 
178a5892487SPeter Brune   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
179ea630c6eSPeter Brune   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
180ea630c6eSPeter Brune   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
181ea630c6eSPeter Brune   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
182ea630c6eSPeter Brune   alpha = alpha / ptAp;
183d9fe6452SJed Brown   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
184a5892487SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
185a5892487SPeter Brune   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
186a5892487SPeter Brune   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
187a5892487SPeter Brune   ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
188ea630c6eSPeter Brune   PetscFunctionReturn(0);
189ea630c6eSPeter Brune }
190ea630c6eSPeter Brune 
19170d3b23bSPeter Brune EXTERN_C_BEGIN
19270d3b23bSPeter Brune #undef __FUNCT__
19370d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NCG"
19470d3b23bSPeter Brune PetscErrorCode  SNESLineSearchSetType_NCG(SNES snes, SNESLineSearchType type)
19570d3b23bSPeter Brune {
19670d3b23bSPeter Brune   PetscErrorCode ierr;
19770d3b23bSPeter Brune   PetscFunctionBegin;
19870d3b23bSPeter Brune 
19970d3b23bSPeter Brune   switch (type) {
20070d3b23bSPeter Brune   case SNES_LS_BASIC:
20170d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
20270d3b23bSPeter Brune     break;
20370d3b23bSPeter Brune   case SNES_LS_BASIC_NONORMS:
20470d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
20570d3b23bSPeter Brune     break;
20670d3b23bSPeter Brune   case SNES_LS_QUADRATIC:
20770d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NCG,PETSC_NULL);CHKERRQ(ierr);
20870d3b23bSPeter Brune     break;
20970d3b23bSPeter Brune   case SNES_LS_TEST:
21070d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchExactLinear_NCG,PETSC_NULL);CHKERRQ(ierr);
21170d3b23bSPeter Brune     break;
212*cf9edfecSPeter Brune   case SNES_LS_SECANT:
213*cf9edfecSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr);
214*cf9edfecSPeter Brune     break;
21570d3b23bSPeter Brune   default:
21670d3b23bSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
21770d3b23bSPeter Brune     break;
21870d3b23bSPeter Brune   }
21970d3b23bSPeter Brune   snes->ls_type = type;
22070d3b23bSPeter Brune   PetscFunctionReturn(0);
22170d3b23bSPeter Brune }
22270d3b23bSPeter Brune EXTERN_C_END
22370d3b23bSPeter Brune 
224fef7b6d8SPeter Brune /*
225fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
226fef7b6d8SPeter Brune 
227fef7b6d8SPeter Brune   Input Parameters:
228fef7b6d8SPeter Brune . snes - the SNES context
229fef7b6d8SPeter Brune 
230fef7b6d8SPeter Brune   Output Parameter:
231fef7b6d8SPeter Brune . outits - number of iterations until termination
232fef7b6d8SPeter Brune 
233fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
234fef7b6d8SPeter Brune */
235fef7b6d8SPeter Brune #undef __FUNCT__
236fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
237fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
238fef7b6d8SPeter Brune {
239a5892487SPeter Brune   Vec                 X, dX, lX, F, W, B, G;
240a5892487SPeter Brune   PetscReal           fnorm, gnorm, dummyNorm;
241fef7b6d8SPeter Brune   PetscScalar         dXdot, dXdot_old;
242fef7b6d8SPeter Brune   PetscInt            maxits, i;
243fef7b6d8SPeter Brune   PetscErrorCode      ierr;
244fef7b6d8SPeter Brune   SNESConvergedReason reason;
245fef7b6d8SPeter Brune   PetscBool           lsSuccess = PETSC_TRUE;
246fef7b6d8SPeter Brune   PetscFunctionBegin;
247fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
248fef7b6d8SPeter Brune 
249fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
250fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
251169e2be8SPeter Brune   dX     = snes->work[1];            /* the steepest direction */
252169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
253fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
254a5892487SPeter Brune   W      = snes->work[0];            /* work vector and solution output for the line search */
255302dbbaeSPeter Brune   B      = snes->vec_rhs;            /* the right hand side */
256a5892487SPeter Brune   G      = snes->work[2];            /* the line search result function evaluation */
257fef7b6d8SPeter Brune 
258fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
259fef7b6d8SPeter Brune   snes->iter = 0;
260fef7b6d8SPeter Brune   snes->norm = 0.;
261fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
262fef7b6d8SPeter Brune 
263fef7b6d8SPeter Brune   /* compute the initial function and preconditioned update delX */
264fef7b6d8SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
265fef7b6d8SPeter Brune   if (snes->domainerror) {
266fef7b6d8SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
267fef7b6d8SPeter Brune     PetscFunctionReturn(0);
268fef7b6d8SPeter Brune   }
269fef7b6d8SPeter Brune   /* convergence test */
270fef7b6d8SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
271fef7b6d8SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
272fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
273fef7b6d8SPeter Brune   snes->norm = fnorm;
274fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
275fef7b6d8SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
276fef7b6d8SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
277fef7b6d8SPeter Brune 
278fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
279fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
280fef7b6d8SPeter Brune   /* test convergence */
281fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
282fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
283fef7b6d8SPeter Brune 
284fef7b6d8SPeter Brune   /* Call general purpose update function */
285fef7b6d8SPeter Brune   if (snes->ops->update) {
286fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
287fef7b6d8SPeter Brune   }
288fef7b6d8SPeter Brune 
289fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
290fef7b6d8SPeter Brune   if (!snes->pc) {
291fef7b6d8SPeter Brune     ierr = VecCopy(F, lX);CHKERRQ(ierr);
292fef7b6d8SPeter Brune     ierr = VecScale(lX, -1.0);CHKERRQ(ierr);
293fef7b6d8SPeter Brune   } else {
294fef7b6d8SPeter Brune     ierr = VecCopy(X, lX);CHKERRQ(ierr);
295cf949ffaSPeter Brune     ierr = SNESSolve(snes->pc, B, lX);CHKERRQ(ierr);
296fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
297302dbbaeSPeter Brune     ierr = VecAXPY(lX, -1.0, X);CHKERRQ(ierr);
29851e86f29SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
299fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
300fef7b6d8SPeter Brune       PetscFunctionReturn(0);
301fef7b6d8SPeter Brune     }
302fef7b6d8SPeter Brune   }
303fef7b6d8SPeter Brune   ierr = VecDot(lX, lX, &dXdot);CHKERRQ(ierr);
304fef7b6d8SPeter Brune   for(i = 1; i < maxits; i++) {
305fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
306a5892487SPeter Brune     ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, lX, fnorm, 0.0, G, W, &dummyNorm, &gnorm, &lsSuccess);CHKERRQ(ierr);
307fef7b6d8SPeter Brune     if (!lsSuccess) {
308fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
309fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
310fef7b6d8SPeter Brune         break;
311fef7b6d8SPeter Brune       }
312fef7b6d8SPeter Brune     }
313fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
314fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
315fef7b6d8SPeter Brune       break;
316fef7b6d8SPeter Brune     }
317fef7b6d8SPeter Brune     if (snes->domainerror) {
318fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
319fef7b6d8SPeter Brune       PetscFunctionReturn(0);
320fef7b6d8SPeter Brune     }
321302dbbaeSPeter Brune 
322a5892487SPeter Brune     /* copy over the solution */
323a5892487SPeter Brune     ierr = VecCopy(G, F);CHKERRQ(ierr);
324a5892487SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
325a5892487SPeter Brune     fnorm = gnorm;
326a5892487SPeter Brune 
327fef7b6d8SPeter Brune     /* Monitor convergence */
328fef7b6d8SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
329169e2be8SPeter Brune     snes->iter = i;
330fef7b6d8SPeter Brune     snes->norm = fnorm;
331fef7b6d8SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
332fef7b6d8SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
333fef7b6d8SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
334302dbbaeSPeter Brune 
335fef7b6d8SPeter Brune     /* Test for convergence */
336fef7b6d8SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
337fef7b6d8SPeter Brune     if (snes->reason) break;
338302dbbaeSPeter Brune 
339302dbbaeSPeter Brune     /* Call general purpose update function */
340302dbbaeSPeter Brune     if (snes->ops->update) {
341302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
342302dbbaeSPeter Brune     }
343ee78dd50SPeter Brune     if (!snes->pc) {
344302dbbaeSPeter Brune       ierr = VecCopy(F,dX);CHKERRQ(ierr);
345302dbbaeSPeter Brune       ierr = VecScale(dX,-1.0);CHKERRQ(ierr);
346302dbbaeSPeter Brune     } else {
347302dbbaeSPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
348cf949ffaSPeter Brune       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
349302dbbaeSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
35051e86f29SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
351302dbbaeSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
352302dbbaeSPeter Brune         PetscFunctionReturn(0);
353302dbbaeSPeter Brune       }
354302dbbaeSPeter Brune       ierr = VecAXPY(dX,-1.0,X);CHKERRQ(ierr);
355302dbbaeSPeter Brune     }
356302dbbaeSPeter Brune 
357169e2be8SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
358302dbbaeSPeter Brune     dXdot_old = dXdot;
359302dbbaeSPeter Brune     ierr = VecDot(dX, dX, &dXdot);CHKERRQ(ierr);
360ee78dd50SPeter Brune #if 0
361421d9b32SPeter Brune     if (1)
362169e2be8SPeter Brune       ierr = PetscPrintf(PETSC_COMM_WORLD, "beta %e = %e / %e \n", dXdot / dXdot_old, dXdot, dXdot_old);CHKERRQ(ierr);
363ee78dd50SPeter Brune #endif
364302dbbaeSPeter Brune     ierr = VecAYPX(lX, dXdot / dXdot_old, dX);CHKERRQ(ierr);
365302dbbaeSPeter Brune     /* line search for the proper contribution of lX to the solution */
366fef7b6d8SPeter Brune   }
367fef7b6d8SPeter Brune   if (i == maxits) {
368fef7b6d8SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
369fef7b6d8SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
370fef7b6d8SPeter Brune   }
371fef7b6d8SPeter Brune   PetscFunctionReturn(0);
372fef7b6d8SPeter Brune }
373fef7b6d8SPeter Brune 
374fef7b6d8SPeter Brune /*MC
375fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
376fef7b6d8SPeter Brune 
377fef7b6d8SPeter Brune   Level: beginner
378fef7b6d8SPeter Brune 
379fef7b6d8SPeter Brune   Options Database:
380fef7b6d8SPeter Brune +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
381fef7b6d8SPeter Brune -   -snes_ls <basic,basicnormnorms,quadratic>
382fef7b6d8SPeter Brune 
383fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
384fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
385fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
386fef7b6d8SPeter Brune 
387fef7b6d8SPeter Brune 
388fef7b6d8SPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
389fef7b6d8SPeter Brune M*/
390fef7b6d8SPeter Brune EXTERN_C_BEGIN
391fef7b6d8SPeter Brune #undef __FUNCT__
392fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
393fef7b6d8SPeter Brune PetscErrorCode  SNESCreate_NCG(SNES snes)
394fef7b6d8SPeter Brune {
395fef7b6d8SPeter Brune   PetscErrorCode   ierr;
396ea630c6eSPeter Brune   SNES_NCG * neP;
397fef7b6d8SPeter Brune 
398fef7b6d8SPeter Brune   PetscFunctionBegin;
399fef7b6d8SPeter Brune   snes->ops->destroy         = SNESDestroy_NCG;
400fef7b6d8SPeter Brune   snes->ops->setup           = SNESSetUp_NCG;
401fef7b6d8SPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
402fef7b6d8SPeter Brune   snes->ops->view            = SNESView_NCG;
403fef7b6d8SPeter Brune   snes->ops->solve           = SNESSolve_NCG;
404fef7b6d8SPeter Brune   snes->ops->reset           = SNESReset_NCG;
405fef7b6d8SPeter Brune 
406fef7b6d8SPeter Brune   snes->usesksp              = PETSC_FALSE;
407fef7b6d8SPeter Brune   snes->usespc               = PETSC_TRUE;
408*cf9edfecSPeter Brune   snes->ls_alpha             = 1e-8;
409fef7b6d8SPeter Brune 
410fef7b6d8SPeter Brune   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
411fef7b6d8SPeter Brune   snes->data = (void*) neP;
41270d3b23bSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NCG",SNESLineSearchSetType_NCG);CHKERRQ(ierr);
413a5892487SPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
414a5892487SPeter Brune 
415fef7b6d8SPeter Brune   PetscFunctionReturn(0);
416fef7b6d8SPeter Brune }
417fef7b6d8SPeter Brune EXTERN_C_END
418