xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 25f8a51bec2b45a4b721a424e44932ad291206ef)
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;
51*dfb256c7SPeter Brune   ierr = SNESDefaultGetWork(snes,4);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 {
66*dfb256c7SPeter Brune   SNES_NCG           *ncg     = (SNES_NCG *)snes->data;
67*dfb256c7SPeter Brune   const char         *betas[] = {"FR","PRP","HS", "DY", "CD"};
68fef7b6d8SPeter Brune   PetscErrorCode     ierr;
69*dfb256c7SPeter Brune   PetscBool          debug, flg;
70*dfb256c7SPeter Brune   PetscInt           indx;
71fef7b6d8SPeter Brune 
72fef7b6d8SPeter Brune   PetscFunctionBegin;
73fef7b6d8SPeter Brune   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
74*dfb256c7SPeter Brune   ierr = PetscOptionsBool("-snes_ncg_monitor", "Monitor NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
75*dfb256c7SPeter Brune   ierr = PetscOptionsEList("-snes_ncg_type","Nonlinear CG update used","", betas, 5, "FR",&indx, &flg);CHKERRQ(ierr);
76*dfb256c7SPeter Brune   if (flg) {
77*dfb256c7SPeter Brune     ncg->betatype = indx;
78*dfb256c7SPeter Brune   }
79*dfb256c7SPeter Brune   if (debug) {
80*dfb256c7SPeter Brune     ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
81*dfb256c7SPeter Brune   }
82fef7b6d8SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
83fef7b6d8SPeter Brune   PetscFunctionReturn(0);
84fef7b6d8SPeter Brune }
85fef7b6d8SPeter Brune 
86fef7b6d8SPeter Brune /*
87fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
88fef7b6d8SPeter Brune 
89fef7b6d8SPeter Brune   Input Parameters:
90fef7b6d8SPeter Brune + SNES - the SNES context
91fef7b6d8SPeter Brune - viewer - visualization context
92fef7b6d8SPeter Brune 
93fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
94fef7b6d8SPeter Brune */
95fef7b6d8SPeter Brune #undef __FUNCT__
96fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
97fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
98fef7b6d8SPeter Brune {
99fef7b6d8SPeter Brune   PetscBool        iascii;
100fef7b6d8SPeter Brune   PetscErrorCode   ierr;
101fef7b6d8SPeter Brune 
102fef7b6d8SPeter Brune   PetscFunctionBegin;
103fef7b6d8SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
104fef7b6d8SPeter Brune   if (iascii) {
105ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
106fef7b6d8SPeter Brune   }
107fef7b6d8SPeter Brune   PetscFunctionReturn(0);
108fef7b6d8SPeter Brune }
109fef7b6d8SPeter Brune 
110fef7b6d8SPeter Brune #undef __FUNCT__
111ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchExactLinear_NCG"
112ea630c6eSPeter 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)
113ea630c6eSPeter Brune {
114ea630c6eSPeter Brune   PetscScalar      alpha, ptAp;
115ea630c6eSPeter Brune   PetscErrorCode   ierr;
116ea630c6eSPeter Brune   /* SNES_NCG *ncg =  (SNES_NCG *) snes->data; */
117ea630c6eSPeter Brune   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
118ea630c6eSPeter Brune 
119ea630c6eSPeter Brune   PetscFunctionBegin;
120ea630c6eSPeter Brune   /*
121ea630c6eSPeter Brune 
122d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
123ea630c6eSPeter Brune 
124ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
125ea630c6eSPeter Brune 
126ea630c6eSPeter Brune    */
127ea630c6eSPeter Brune 
128a5892487SPeter Brune   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
129ea630c6eSPeter Brune   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
130ea630c6eSPeter Brune   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
131ea630c6eSPeter Brune   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
132ea630c6eSPeter Brune   alpha = alpha / ptAp;
133d9fe6452SJed Brown   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
134a5892487SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
135a5892487SPeter Brune   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
136a5892487SPeter Brune   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
137a5892487SPeter Brune   ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
138ea630c6eSPeter Brune   PetscFunctionReturn(0);
139ea630c6eSPeter Brune }
140ea630c6eSPeter Brune 
14170d3b23bSPeter Brune EXTERN_C_BEGIN
14270d3b23bSPeter Brune #undef __FUNCT__
14370d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NCG"
14470d3b23bSPeter Brune PetscErrorCode  SNESLineSearchSetType_NCG(SNES snes, SNESLineSearchType type)
14570d3b23bSPeter Brune {
14670d3b23bSPeter Brune   PetscErrorCode ierr;
14770d3b23bSPeter Brune   PetscFunctionBegin;
14870d3b23bSPeter Brune 
14970d3b23bSPeter Brune   switch (type) {
15070d3b23bSPeter Brune   case SNES_LS_BASIC:
15170d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
15270d3b23bSPeter Brune     break;
15370d3b23bSPeter Brune   case SNES_LS_BASIC_NONORMS:
15470d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
15570d3b23bSPeter Brune     break;
15670d3b23bSPeter Brune   case SNES_LS_QUADRATIC:
157d2140ca3SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
15870d3b23bSPeter Brune     break;
15970d3b23bSPeter Brune   case SNES_LS_TEST:
16070d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchExactLinear_NCG,PETSC_NULL);CHKERRQ(ierr);
16170d3b23bSPeter Brune     break;
162cf9edfecSPeter Brune   case SNES_LS_SECANT:
163cf9edfecSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr);
164cf9edfecSPeter Brune     break;
16570d3b23bSPeter Brune   default:
16670d3b23bSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
16770d3b23bSPeter Brune     break;
16870d3b23bSPeter Brune   }
16970d3b23bSPeter Brune   snes->ls_type = type;
17070d3b23bSPeter Brune   PetscFunctionReturn(0);
17170d3b23bSPeter Brune }
17270d3b23bSPeter Brune EXTERN_C_END
17370d3b23bSPeter Brune 
174fef7b6d8SPeter Brune /*
175fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
176fef7b6d8SPeter Brune 
177fef7b6d8SPeter Brune   Input Parameters:
178fef7b6d8SPeter Brune . snes - the SNES context
179fef7b6d8SPeter Brune 
180fef7b6d8SPeter Brune   Output Parameter:
181fef7b6d8SPeter Brune . outits - number of iterations until termination
182fef7b6d8SPeter Brune 
183fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
184fef7b6d8SPeter Brune */
185fef7b6d8SPeter Brune #undef __FUNCT__
186fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
187fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
188fef7b6d8SPeter Brune {
189*dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
19029c94e34SPeter Brune   Vec                 X, dX, lX, F, W, B, G, dXold;
191*dfb256c7SPeter Brune   PetscReal           fnorm, gnorm, dummyNorm, beta = 0.0;
19229c94e34SPeter Brune   PetscScalar         dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
193fef7b6d8SPeter Brune   PetscInt            maxits, i;
194fef7b6d8SPeter Brune   PetscErrorCode      ierr;
195fef7b6d8SPeter Brune   SNESConvergedReason reason;
196fef7b6d8SPeter Brune   PetscBool           lsSuccess = PETSC_TRUE;
197fef7b6d8SPeter Brune   PetscFunctionBegin;
198fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
199fef7b6d8SPeter Brune 
200fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
201fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
20229c94e34SPeter Brune   dXold  = snes->work[0];            /* The previous iterate of X */
203169e2be8SPeter Brune   dX     = snes->work[1];            /* the steepest direction */
204169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
205fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
206*dfb256c7SPeter Brune   W      = snes->work[2];            /* work vector and solution output for the line search */
207302dbbaeSPeter Brune   B      = snes->vec_rhs;            /* the right hand side */
208*dfb256c7SPeter Brune   G      = snes->work[3];            /* the line search result function evaluation */
209fef7b6d8SPeter Brune 
210fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
211fef7b6d8SPeter Brune   snes->iter = 0;
212fef7b6d8SPeter Brune   snes->norm = 0.;
213fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
214fef7b6d8SPeter Brune 
215fef7b6d8SPeter Brune   /* compute the initial function and preconditioned update delX */
216fef7b6d8SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
217fef7b6d8SPeter Brune   if (snes->domainerror) {
218fef7b6d8SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
219fef7b6d8SPeter Brune     PetscFunctionReturn(0);
220fef7b6d8SPeter Brune   }
221fef7b6d8SPeter Brune   /* convergence test */
222fef7b6d8SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
223fef7b6d8SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
224fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
225fef7b6d8SPeter Brune   snes->norm = fnorm;
226fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
227fef7b6d8SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
228fef7b6d8SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
229fef7b6d8SPeter Brune 
230fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
231fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
232fef7b6d8SPeter Brune   /* test convergence */
233fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
234fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
235fef7b6d8SPeter Brune 
236fef7b6d8SPeter Brune   /* Call general purpose update function */
237fef7b6d8SPeter Brune   if (snes->ops->update) {
238fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
239fef7b6d8SPeter Brune  }
240fef7b6d8SPeter Brune 
241fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
2425d10c001SPeter Brune   if (snes->usegs && snes->ops->computegs) {
24329c94e34SPeter Brune     ierr = VecCopy(X, dX);CHKERRQ(ierr);
24429c94e34SPeter Brune     ierr = SNESComputeGS(snes, B, dX);CHKERRQ(ierr);
24529c94e34SPeter Brune     ierr = VecAYPX(dX, -1.0, X);CHKERRQ(ierr);
2465d10c001SPeter Brune   } else if (snes->pc) {
24729c94e34SPeter Brune     ierr = VecCopy(X, dX);CHKERRQ(ierr);
24829c94e34SPeter Brune     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
249fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
25051e86f29SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
251fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
252fef7b6d8SPeter Brune       PetscFunctionReturn(0);
253fef7b6d8SPeter Brune     }
25429c94e34SPeter Brune     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
2555d10c001SPeter Brune   } else {
25629c94e34SPeter Brune     ierr = VecCopy(F, dX);CHKERRQ(ierr);
257fef7b6d8SPeter Brune   }
25829c94e34SPeter Brune   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
25929c94e34SPeter Brune   ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
26009c08436SPeter Brune   for(i = 1; i < maxits + 1; i++) {
261fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
26229c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
26329c94e34SPeter Brune     if (ncg->betatype == 1 || ncg->betatype == 2 || ncg->betatype == 3 || ncg->betatype == 4) { /* prp, hs, dy, cd need dXold*/
26429c94e34SPeter Brune       ierr = VecCopy(dX, dXold);CHKERRQ(ierr);
26529c94e34SPeter Brune     }
266a5892487SPeter Brune     ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, lX, fnorm, 0.0, G, W, &dummyNorm, &gnorm, &lsSuccess);CHKERRQ(ierr);
267fef7b6d8SPeter Brune     if (!lsSuccess) {
268fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
269fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2705d10c001SPeter Brune         PetscFunctionReturn(0);
271fef7b6d8SPeter Brune       }
272fef7b6d8SPeter Brune     }
273fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
274fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2755d10c001SPeter Brune       PetscFunctionReturn(0);
276fef7b6d8SPeter Brune     }
277fef7b6d8SPeter Brune     if (snes->domainerror) {
278fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
279fef7b6d8SPeter Brune       PetscFunctionReturn(0);
280fef7b6d8SPeter Brune     }
281302dbbaeSPeter Brune 
282a5892487SPeter Brune     /* copy over the solution */
283a5892487SPeter Brune     ierr = VecCopy(G, F);CHKERRQ(ierr);
284a5892487SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
285a5892487SPeter Brune     fnorm = gnorm;
286a5892487SPeter Brune 
287fef7b6d8SPeter Brune     /* Monitor convergence */
288fef7b6d8SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
289169e2be8SPeter Brune     snes->iter = i;
290fef7b6d8SPeter Brune     snes->norm = fnorm;
291fef7b6d8SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
292fef7b6d8SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
293fef7b6d8SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
294302dbbaeSPeter Brune 
295fef7b6d8SPeter Brune     /* Test for convergence */
296fef7b6d8SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2975d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
298302dbbaeSPeter Brune 
299302dbbaeSPeter Brune     /* Call general purpose update function */
300302dbbaeSPeter Brune     if (snes->ops->update) {
301302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
302302dbbaeSPeter Brune     }
303a3685310SPeter Brune     if (snes->usegs && snes->ops->computegs) {
304a3685310SPeter Brune       ierr = VecCopy(X, dX);CHKERRQ(ierr);
305a3685310SPeter Brune       ierr = SNESComputeGS(snes, B, dX);CHKERRQ(ierr);
30629c94e34SPeter Brune       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
307a3685310SPeter Brune     } else if (snes->pc) {
308302dbbaeSPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
309cf949ffaSPeter Brune       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
310302dbbaeSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
31151e86f29SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
312302dbbaeSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
313302dbbaeSPeter Brune         PetscFunctionReturn(0);
314302dbbaeSPeter Brune       }
315eb1825c3SPeter Brune       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
316a3685310SPeter Brune     } else {
317a3685310SPeter Brune       ierr = VecCopy(F,dX);CHKERRQ(ierr);
318302dbbaeSPeter Brune     }
319302dbbaeSPeter Brune 
32029c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
321*dfb256c7SPeter Brune     switch(ncg->betatype) {
322*dfb256c7SPeter Brune     case 0: /* Fletcher-Reeves */
32329c94e34SPeter Brune       dXolddotdXold = dXdotdX;
32429c94e34SPeter Brune       ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
32529c94e34SPeter Brune       beta = PetscRealPart(dXdotdX / dXolddotdXold);
326*dfb256c7SPeter Brune       break;
327*dfb256c7SPeter Brune     case 1: /* Polak-Ribiere-Poylak */
32829c94e34SPeter Brune       dXolddotdXold = dXdotdX;
32929c94e34SPeter Brune       ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
33029c94e34SPeter Brune       ierr = VecDot(dXold, dX, &dXdotdXold);CHKERRQ(ierr);
33125f8a51bSPeter Brune       beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
332*dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
333*dfb256c7SPeter Brune       break;
334*dfb256c7SPeter Brune     case 2: /* Hestenes-Stiefel */
33529c94e34SPeter Brune       ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
33629c94e34SPeter Brune       ierr = VecDot(dX, dXold, &dXdotdXold);CHKERRQ(ierr);
33729c94e34SPeter Brune       ierr = VecDot(lX, dX, &lXdotdX);CHKERRQ(ierr);
33829c94e34SPeter Brune       ierr = VecDot(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
33929c94e34SPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
340*dfb256c7SPeter Brune       break;
341*dfb256c7SPeter Brune     case 3: /* Dai-Yuan */
34229c94e34SPeter Brune       ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
34329c94e34SPeter Brune       ierr = VecDot(lX, dX, &lXdotdXold);CHKERRQ(ierr);
34429c94e34SPeter Brune       ierr = VecDot(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
34529c94e34SPeter Brune       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr);
346*dfb256c7SPeter Brune       break;
347*dfb256c7SPeter Brune     case 4: /* Conjugate Descent */
34829c94e34SPeter Brune       ierr = VecDot(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
34929c94e34SPeter Brune       ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
35029c94e34SPeter Brune       beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr);
351*dfb256c7SPeter Brune       break;
352*dfb256c7SPeter Brune     }
353*dfb256c7SPeter Brune     if (ncg->monitor) {
354646217ecSPeter Brune       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
355*dfb256c7SPeter Brune     }
356*dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
357fef7b6d8SPeter Brune   }
358fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
359fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
360fef7b6d8SPeter Brune   PetscFunctionReturn(0);
361fef7b6d8SPeter Brune }
362fef7b6d8SPeter Brune 
363fef7b6d8SPeter Brune /*MC
364fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
365fef7b6d8SPeter Brune 
366fef7b6d8SPeter Brune   Level: beginner
367fef7b6d8SPeter Brune 
368fef7b6d8SPeter Brune   Options Database:
369fef7b6d8SPeter Brune +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
370fef7b6d8SPeter Brune -   -snes_ls <basic,basicnormnorms,quadratic>
371fef7b6d8SPeter Brune 
372fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
373fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
374fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
375fef7b6d8SPeter Brune 
376fef7b6d8SPeter Brune 
377fef7b6d8SPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
378fef7b6d8SPeter Brune M*/
379fef7b6d8SPeter Brune EXTERN_C_BEGIN
380fef7b6d8SPeter Brune #undef __FUNCT__
381fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
382fef7b6d8SPeter Brune PetscErrorCode  SNESCreate_NCG(SNES snes)
383fef7b6d8SPeter Brune {
384fef7b6d8SPeter Brune   PetscErrorCode   ierr;
385ea630c6eSPeter Brune   SNES_NCG * neP;
386fef7b6d8SPeter Brune 
387fef7b6d8SPeter Brune   PetscFunctionBegin;
388fef7b6d8SPeter Brune   snes->ops->destroy         = SNESDestroy_NCG;
389fef7b6d8SPeter Brune   snes->ops->setup           = SNESSetUp_NCG;
390fef7b6d8SPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
391fef7b6d8SPeter Brune   snes->ops->view            = SNESView_NCG;
392fef7b6d8SPeter Brune   snes->ops->solve           = SNESSolve_NCG;
393fef7b6d8SPeter Brune   snes->ops->reset           = SNESReset_NCG;
394fef7b6d8SPeter Brune 
395fef7b6d8SPeter Brune   snes->usesksp              = PETSC_FALSE;
396fef7b6d8SPeter Brune   snes->usespc               = PETSC_TRUE;
397fef7b6d8SPeter Brune 
398fef7b6d8SPeter Brune   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
399fef7b6d8SPeter Brune   snes->data = (void*) neP;
400*dfb256c7SPeter Brune   neP->monitor = PETSC_NULL;
40125f8a51bSPeter Brune   neP->betatype = 4;
40270d3b23bSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NCG",SNESLineSearchSetType_NCG);CHKERRQ(ierr);
403d2140ca3SPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
404fef7b6d8SPeter Brune   PetscFunctionReturn(0);
405fef7b6d8SPeter Brune }
406fef7b6d8SPeter Brune EXTERN_C_END
407