xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 1867fe5b73a024d1fa5fcd64303088590170b736)
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 
8fef7b6d8SPeter Brune   PetscFunctionBegin;
9fef7b6d8SPeter Brune   PetscFunctionReturn(0);
10fef7b6d8SPeter Brune }
11fef7b6d8SPeter Brune 
12fef7b6d8SPeter Brune /*
13fef7b6d8SPeter Brune   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
14fef7b6d8SPeter Brune 
15fef7b6d8SPeter Brune   Input Parameter:
16fef7b6d8SPeter Brune . snes - the SNES context
17fef7b6d8SPeter Brune 
18fef7b6d8SPeter Brune   Application Interface Routine: SNESDestroy()
19fef7b6d8SPeter Brune */
20fef7b6d8SPeter Brune #undef __FUNCT__
21fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG"
22fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes)
23fef7b6d8SPeter Brune {
24fef7b6d8SPeter Brune   PetscErrorCode   ierr;
25fef7b6d8SPeter Brune 
26fef7b6d8SPeter Brune   PetscFunctionBegin;
27fef7b6d8SPeter Brune   ierr = SNESReset_NCG(snes);CHKERRQ(ierr);
28fef7b6d8SPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
29fef7b6d8SPeter Brune   PetscFunctionReturn(0);
30fef7b6d8SPeter Brune }
31fef7b6d8SPeter Brune 
32fef7b6d8SPeter Brune /*
33fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
34fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
35fef7b6d8SPeter Brune 
36fef7b6d8SPeter Brune    Input Parameters:
37fef7b6d8SPeter Brune +  snes - the SNES context
38fef7b6d8SPeter Brune -  x - the solution vector
39fef7b6d8SPeter Brune 
40fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
41fef7b6d8SPeter Brune  */
42fef7b6d8SPeter Brune #undef __FUNCT__
43fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG"
44fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes)
45fef7b6d8SPeter Brune {
46fef7b6d8SPeter Brune   PetscErrorCode ierr;
47fef7b6d8SPeter Brune 
48fef7b6d8SPeter Brune   PetscFunctionBegin;
4988d374b2SPeter Brune   ierr = SNESDefaultGetWork(snes,5);CHKERRQ(ierr);
50fef7b6d8SPeter Brune   PetscFunctionReturn(0);
51fef7b6d8SPeter Brune }
52fef7b6d8SPeter Brune /*
5305b53524SPeter Brune   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.
54fef7b6d8SPeter Brune 
55fef7b6d8SPeter Brune   Input Parameter:
56fef7b6d8SPeter Brune . snes - the SNES context
57fef7b6d8SPeter Brune 
58fef7b6d8SPeter Brune   Application Interface Routine: SNESSetFromOptions()
59fef7b6d8SPeter Brune */
60fef7b6d8SPeter Brune #undef __FUNCT__
61fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG"
62fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
63fef7b6d8SPeter Brune {
64dfb256c7SPeter Brune   SNES_NCG           *ncg     = (SNES_NCG *)snes->data;
65dfb256c7SPeter Brune   const char         *betas[] = {"FR","PRP","HS", "DY", "CD"};
66fef7b6d8SPeter Brune   PetscErrorCode     ierr;
67dfb256c7SPeter Brune   PetscBool          debug, flg;
68dfb256c7SPeter Brune   PetscInt           indx;
69fef7b6d8SPeter Brune 
70fef7b6d8SPeter Brune   PetscFunctionBegin;
71fef7b6d8SPeter Brune   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
72dfb256c7SPeter Brune   ierr = PetscOptionsBool("-snes_ncg_monitor", "Monitor NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr);
73dfb256c7SPeter Brune   ierr = PetscOptionsEList("-snes_ncg_type","Nonlinear CG update used","", betas, 5, "FR",&indx, &flg);CHKERRQ(ierr);
74dfb256c7SPeter Brune   if (flg) {
75dfb256c7SPeter Brune     ncg->betatype = indx;
76dfb256c7SPeter Brune   }
77dfb256c7SPeter Brune   if (debug) {
78dfb256c7SPeter Brune     ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
79dfb256c7SPeter Brune   }
80fef7b6d8SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
81fef7b6d8SPeter Brune   PetscFunctionReturn(0);
82fef7b6d8SPeter Brune }
83fef7b6d8SPeter Brune 
84fef7b6d8SPeter Brune /*
85fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
86fef7b6d8SPeter Brune 
87fef7b6d8SPeter Brune   Input Parameters:
88fef7b6d8SPeter Brune + SNES - the SNES context
89fef7b6d8SPeter Brune - viewer - visualization context
90fef7b6d8SPeter Brune 
91fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
92fef7b6d8SPeter Brune */
93fef7b6d8SPeter Brune #undef __FUNCT__
94fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
95fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
96fef7b6d8SPeter Brune {
97fef7b6d8SPeter Brune   PetscBool        iascii;
98fef7b6d8SPeter Brune   PetscErrorCode   ierr;
99fef7b6d8SPeter Brune 
100fef7b6d8SPeter Brune   PetscFunctionBegin;
101fef7b6d8SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
102fef7b6d8SPeter Brune   if (iascii) {
103ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
104fef7b6d8SPeter Brune   }
105fef7b6d8SPeter Brune   PetscFunctionReturn(0);
106fef7b6d8SPeter Brune }
107fef7b6d8SPeter Brune 
108fef7b6d8SPeter Brune #undef __FUNCT__
109ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchExactLinear_NCG"
110ea630c6eSPeter 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)
111ea630c6eSPeter Brune {
112ea630c6eSPeter Brune   PetscScalar      alpha, ptAp;
113ea630c6eSPeter Brune   PetscErrorCode   ierr;
114ea630c6eSPeter Brune   /* SNES_NCG *ncg =  (SNES_NCG *) snes->data; */
115ea630c6eSPeter Brune   MatStructure     flg = DIFFERENT_NONZERO_PATTERN;
116ea630c6eSPeter Brune 
117ea630c6eSPeter Brune   PetscFunctionBegin;
118ea630c6eSPeter Brune   /*
119ea630c6eSPeter Brune 
120d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
121ea630c6eSPeter Brune 
122ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
123ea630c6eSPeter Brune 
124ea630c6eSPeter Brune    */
125ea630c6eSPeter Brune 
126a5892487SPeter Brune   ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
127ea630c6eSPeter Brune   ierr = VecDot(F, F, &alpha);CHKERRQ(ierr);
128ea630c6eSPeter Brune   ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
129ea630c6eSPeter Brune   ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
130ea630c6eSPeter Brune   alpha = alpha / ptAp;
131d9fe6452SJed Brown   ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
132a5892487SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
133a5892487SPeter Brune   ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr);
134a5892487SPeter Brune   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
135a5892487SPeter Brune   ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
136ea630c6eSPeter Brune   PetscFunctionReturn(0);
137ea630c6eSPeter Brune }
138ea630c6eSPeter Brune 
13988d374b2SPeter Brune #undef __FUNCT__
14088d374b2SPeter Brune #define __FUNCT__ "ComputeYtJtF_Private"
14188d374b2SPeter Brune /*
14288d374b2SPeter Brune 
14388d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
14488d374b2SPeter Brune 
14588d374b2SPeter Brune  */
14688d374b2SPeter Brune PetscErrorCode ComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) {
14788d374b2SPeter Brune   PetscErrorCode ierr;
14888d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
14988d374b2SPeter Brune   PetscFunctionBegin;
15088d374b2SPeter Brune   ierr = VecDot(F, F, &ftf);CHKERRQ(ierr);
15188d374b2SPeter Brune   ierr = VecDot(F, Y, &fty);CHKERRQ(ierr);
15288d374b2SPeter Brune   h = 1e-5*fty / fty;
15388d374b2SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
15488d374b2SPeter Brune   ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr);            /* this is arbitrary */
15588d374b2SPeter Brune   ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
15688d374b2SPeter Brune   ierr = VecDot(G, F, &ftg);CHKERRQ(ierr);
15788d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
15888d374b2SPeter Brune   PetscFunctionReturn(0);
15988d374b2SPeter Brune }
16088d374b2SPeter Brune 
16170d3b23bSPeter Brune EXTERN_C_BEGIN
16270d3b23bSPeter Brune #undef __FUNCT__
16370d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NCG"
16470d3b23bSPeter Brune PetscErrorCode  SNESLineSearchSetType_NCG(SNES snes, SNESLineSearchType type)
16570d3b23bSPeter Brune {
16670d3b23bSPeter Brune   PetscErrorCode ierr;
16770d3b23bSPeter Brune   PetscFunctionBegin;
16870d3b23bSPeter Brune 
16970d3b23bSPeter Brune   switch (type) {
17070d3b23bSPeter Brune   case SNES_LS_BASIC:
17170d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
17270d3b23bSPeter Brune     break;
17370d3b23bSPeter Brune   case SNES_LS_BASIC_NONORMS:
17470d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
17570d3b23bSPeter Brune     break;
17670d3b23bSPeter Brune   case SNES_LS_QUADRATIC:
177d2140ca3SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
17870d3b23bSPeter Brune     break;
17970d3b23bSPeter Brune   case SNES_LS_TEST:
18070d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchExactLinear_NCG,PETSC_NULL);CHKERRQ(ierr);
18170d3b23bSPeter Brune     break;
182cf9edfecSPeter Brune   case SNES_LS_SECANT:
183cf9edfecSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr);
184cf9edfecSPeter Brune     break;
18570d3b23bSPeter Brune   default:
18670d3b23bSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
18770d3b23bSPeter Brune     break;
18870d3b23bSPeter Brune   }
18970d3b23bSPeter Brune   snes->ls_type = type;
19070d3b23bSPeter Brune   PetscFunctionReturn(0);
19170d3b23bSPeter Brune }
19270d3b23bSPeter Brune EXTERN_C_END
19370d3b23bSPeter Brune 
194fef7b6d8SPeter Brune /*
195fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
196fef7b6d8SPeter Brune 
197fef7b6d8SPeter Brune   Input Parameters:
198fef7b6d8SPeter Brune . snes - the SNES context
199fef7b6d8SPeter Brune 
200fef7b6d8SPeter Brune   Output Parameter:
201fef7b6d8SPeter Brune . outits - number of iterations until termination
202fef7b6d8SPeter Brune 
203fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
204fef7b6d8SPeter Brune */
205fef7b6d8SPeter Brune #undef __FUNCT__
206fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
207fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
208fef7b6d8SPeter Brune {
209dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
21088d374b2SPeter Brune   Vec                 X, dX, lX, F, W, B, G, Fold, Xold;
211dfb256c7SPeter Brune   PetscReal           fnorm, gnorm, dummyNorm, beta = 0.0;
2125d115551SPeter Brune   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
213fef7b6d8SPeter Brune   PetscInt            maxits, i;
214fef7b6d8SPeter Brune   PetscErrorCode      ierr;
215fef7b6d8SPeter Brune   SNESConvergedReason reason;
216fef7b6d8SPeter Brune   PetscBool           lsSuccess = PETSC_TRUE;
21788d374b2SPeter Brune 
218fef7b6d8SPeter Brune   PetscFunctionBegin;
219fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
220fef7b6d8SPeter Brune 
221fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
222fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
2235d115551SPeter Brune   Fold   = snes->work[0];            /* The previous iterate of X */
22488d374b2SPeter Brune   dX     = snes->work[1];            /* the preconditioned direction */
225169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
226fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
227dfb256c7SPeter Brune   W      = snes->work[2];            /* work vector and solution output for the line search */
228302dbbaeSPeter Brune   B      = snes->vec_rhs;            /* the right hand side */
229dfb256c7SPeter Brune   G      = snes->work[3];            /* the line search result function evaluation */
23088d374b2SPeter Brune   Xold   = snes->work[4];            /* the last solution (necessary for quadratic line search to not stagnate) */
231fef7b6d8SPeter Brune 
232fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
233fef7b6d8SPeter Brune   snes->iter = 0;
234fef7b6d8SPeter Brune   snes->norm = 0.;
235fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
236fef7b6d8SPeter Brune 
237fef7b6d8SPeter Brune   /* compute the initial function and preconditioned update delX */
238fef7b6d8SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
239fef7b6d8SPeter Brune   if (snes->domainerror) {
240fef7b6d8SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
241fef7b6d8SPeter Brune     PetscFunctionReturn(0);
242fef7b6d8SPeter Brune   }
243fef7b6d8SPeter Brune   /* convergence test */
244fef7b6d8SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
245fef7b6d8SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
246fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
247fef7b6d8SPeter Brune   snes->norm = fnorm;
248fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
249fef7b6d8SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
250fef7b6d8SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
251fef7b6d8SPeter Brune 
252fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
253fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
254fef7b6d8SPeter Brune   /* test convergence */
255fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
256fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
257fef7b6d8SPeter Brune 
258fef7b6d8SPeter Brune   /* Call general purpose update function */
259fef7b6d8SPeter Brune   if (snes->ops->update) {
260fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
261fef7b6d8SPeter Brune  }
262fef7b6d8SPeter Brune 
263fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
26488d374b2SPeter Brune   /*
26588d374b2SPeter Brune   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
26688d374b2SPeter Brune   ierr = MatMultTranspose(snes->jacobian, F, dF);CHKERRQ(ierr);
26788d374b2SPeter Brune    */
26883947a81SPeter Brune   if (snes->pc) {
26929c94e34SPeter Brune     ierr = VecCopy(X, dX);CHKERRQ(ierr);
27029c94e34SPeter Brune     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
271fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
27251e86f29SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
273fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
274fef7b6d8SPeter Brune       PetscFunctionReturn(0);
275fef7b6d8SPeter Brune     }
27629c94e34SPeter Brune     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
2775d10c001SPeter Brune   } else {
27829c94e34SPeter Brune     ierr = VecCopy(F, dX);CHKERRQ(ierr);
279fef7b6d8SPeter Brune   }
28088d374b2SPeter Brune 
28129c94e34SPeter Brune   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
28288d374b2SPeter Brune 
28388d374b2SPeter Brune   if (snes->ls_type == SNES_LS_SECANT) {
2845d115551SPeter Brune     ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
28588d374b2SPeter Brune   } else {
28688d374b2SPeter Brune     ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
28788d374b2SPeter Brune   }
2885d115551SPeter Brune 
28909c08436SPeter Brune   for(i = 1; i < maxits + 1; i++) {
290fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
29129c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
29229c94e34SPeter Brune     if (ncg->betatype == 1 || ncg->betatype == 2 || ncg->betatype == 3 || ncg->betatype == 4) { /* prp, hs, dy, cd need dXold*/
2935d115551SPeter Brune       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
29488d374b2SPeter Brune       if (snes->ls_type == SNES_LS_QUADRATIC) {
29588d374b2SPeter Brune         ierr = VecCopy(X, Xold);CHKERRQ(ierr);
29688d374b2SPeter Brune       }
29729c94e34SPeter Brune     }
29805b53524SPeter Brune     ierr = SNESLineSearchPreCheckApply(snes, X, lX, PETSC_NULL);CHKERRQ(ierr);
29905b53524SPeter Brune     ierr = SNESLineSearchApply(snes, X, F, lX, fnorm, 0.0, W, G, &dummyNorm, &gnorm, &lsSuccess);CHKERRQ(ierr);
300fef7b6d8SPeter Brune     if (!lsSuccess) {
301fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
302fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3035d10c001SPeter Brune         PetscFunctionReturn(0);
304fef7b6d8SPeter Brune       }
305fef7b6d8SPeter Brune     }
306fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
307fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3085d10c001SPeter Brune       PetscFunctionReturn(0);
309fef7b6d8SPeter Brune     }
310fef7b6d8SPeter Brune     if (snes->domainerror) {
311fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
312fef7b6d8SPeter Brune       PetscFunctionReturn(0);
313fef7b6d8SPeter Brune     }
314302dbbaeSPeter Brune 
315a5892487SPeter Brune     /* copy over the solution */
316a5892487SPeter Brune     ierr = VecCopy(G, F);CHKERRQ(ierr);
317a5892487SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
318a5892487SPeter Brune     fnorm = gnorm;
319a5892487SPeter Brune 
320fef7b6d8SPeter Brune     /* Monitor convergence */
321fef7b6d8SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
322169e2be8SPeter Brune     snes->iter = i;
323fef7b6d8SPeter Brune     snes->norm = fnorm;
324fef7b6d8SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
325fef7b6d8SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
326fef7b6d8SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
327302dbbaeSPeter Brune 
328fef7b6d8SPeter Brune     /* Test for convergence */
329fef7b6d8SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3305d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
331302dbbaeSPeter Brune 
332302dbbaeSPeter Brune     /* Call general purpose update function */
333302dbbaeSPeter Brune     if (snes->ops->update) {
334302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
335302dbbaeSPeter Brune     }
33683947a81SPeter Brune     if (snes->pc) {
337302dbbaeSPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
338cf949ffaSPeter Brune       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
339302dbbaeSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
34051e86f29SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
341302dbbaeSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
342302dbbaeSPeter Brune         PetscFunctionReturn(0);
343302dbbaeSPeter Brune       }
344eb1825c3SPeter Brune       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
345a3685310SPeter Brune     } else {
346a3685310SPeter Brune       ierr = VecCopy(F, dX);CHKERRQ(ierr);
347302dbbaeSPeter Brune     }
348302dbbaeSPeter Brune 
34929c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
350dfb256c7SPeter Brune     switch(ncg->betatype) {
351dfb256c7SPeter Brune     case 0: /* Fletcher-Reeves */
3525d115551SPeter Brune       dXolddotFold = dXdotF;
35388d374b2SPeter Brune       if (snes->ls_type == SNES_LS_SECANT) {
3545d115551SPeter Brune         ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
35588d374b2SPeter Brune       } else {
35688d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
3575d115551SPeter Brune         beta = PetscRealPart(dXdotF / dXolddotFold);
35888d374b2SPeter Brune       }
359dfb256c7SPeter Brune       break;
360dfb256c7SPeter Brune     case 1: /* Polak-Ribiere-Poylak */
3615d115551SPeter Brune       dXolddotFold = dXdotF;
36288d374b2SPeter Brune       if (snes->ls_type == SNES_LS_SECANT) {
3635d115551SPeter Brune         ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
3645d115551SPeter Brune         ierr = VecDot(Fold, dX, &dXdotFold);CHKERRQ(ierr);
36588d374b2SPeter Brune       } else {
36688d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
36788d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, Xold, Fold, dX, W, G, &dXdotFold);CHKERRQ(ierr);
36888d374b2SPeter Brune       }
3695d115551SPeter Brune       beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
370dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
371dfb256c7SPeter Brune       break;
372dfb256c7SPeter Brune     case 2: /* Hestenes-Stiefel */
37388d374b2SPeter Brune       if (snes->ls_type == SNES_LS_SECANT) {
3745d115551SPeter Brune         ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
3755d115551SPeter Brune         ierr = VecDot(dX, Fold, &dXdotFold);CHKERRQ(ierr);
3765d115551SPeter Brune         ierr = VecDot(lX, F, &lXdotF);CHKERRQ(ierr);
3775d115551SPeter Brune         ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
37888d374b2SPeter Brune       } else {
37988d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
38088d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, Xold, Fold, dX, W, G, &dXdotFold);CHKERRQ(ierr);
38188d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, lX, W, G, &lXdotF);CHKERRQ(ierr);
38288d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, Xold, Fold, lX, W, G, &lXdotFold);CHKERRQ(ierr);
38388d374b2SPeter Brune       }
3845d115551SPeter Brune       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
385dfb256c7SPeter Brune       break;
386dfb256c7SPeter Brune     case 3: /* Dai-Yuan */
38788d374b2SPeter Brune       if (snes->ls_type == SNES_LS_SECANT) {
3885d115551SPeter Brune         ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
38988d374b2SPeter Brune         ierr = VecDot(dX, F, &lXdotF);CHKERRQ(ierr);
39088d374b2SPeter Brune         ierr = VecDot(dX, F, &lXdotFold);CHKERRQ(ierr);
39188d374b2SPeter Brune       } else {
39288d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
39388d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, lX, W, G, &dXdotF);CHKERRQ(ierr);
39488d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, Xold, Fold, lX, W, G, &lXdotFold);CHKERRQ(ierr);
39588d374b2SPeter Brune       }
3965d115551SPeter Brune       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
397dfb256c7SPeter Brune       break;
398dfb256c7SPeter Brune     case 4: /* Conjugate Descent */
39988d374b2SPeter Brune       if (snes->ls_type == SNES_LS_SECANT) {
40001fcfa61SPeter Brune         ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr);
40188d374b2SPeter Brune         ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr);
40288d374b2SPeter Brune       } else {
40388d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
40488d374b2SPeter Brune         ierr = ComputeYtJtF_Private(snes, Xold, Fold, lX, W, G, &lXdotFold);CHKERRQ(ierr);
40588d374b2SPeter Brune       }
4065d115551SPeter Brune       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
407dfb256c7SPeter Brune       break;
408dfb256c7SPeter Brune     }
409dfb256c7SPeter Brune     if (ncg->monitor) {
410646217ecSPeter Brune       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
411dfb256c7SPeter Brune     }
412dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
413fef7b6d8SPeter Brune   }
414fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
415fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
416fef7b6d8SPeter Brune   PetscFunctionReturn(0);
417fef7b6d8SPeter Brune }
418fef7b6d8SPeter Brune 
419fef7b6d8SPeter Brune /*MC
420fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
421fef7b6d8SPeter Brune 
422fef7b6d8SPeter Brune   Level: beginner
423fef7b6d8SPeter Brune 
424fef7b6d8SPeter Brune   Options Database:
425*1867fe5bSPeter Brune +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
426*1867fe5bSPeter Brune .   -snes_ls <basic,basicnormnorms,quadratic,secant,test> - Line search type.
427*1867fe5bSPeter Brune -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
428fef7b6d8SPeter Brune 
429fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
430fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
431fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
432fef7b6d8SPeter Brune 
433fef7b6d8SPeter Brune 
434fef7b6d8SPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
435fef7b6d8SPeter Brune M*/
436fef7b6d8SPeter Brune EXTERN_C_BEGIN
437fef7b6d8SPeter Brune #undef __FUNCT__
438fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
439fef7b6d8SPeter Brune PetscErrorCode  SNESCreate_NCG(SNES snes)
440fef7b6d8SPeter Brune {
441fef7b6d8SPeter Brune   PetscErrorCode   ierr;
442ea630c6eSPeter Brune   SNES_NCG * neP;
443fef7b6d8SPeter Brune 
444fef7b6d8SPeter Brune   PetscFunctionBegin;
445fef7b6d8SPeter Brune   snes->ops->destroy         = SNESDestroy_NCG;
446fef7b6d8SPeter Brune   snes->ops->setup           = SNESSetUp_NCG;
447fef7b6d8SPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
448fef7b6d8SPeter Brune   snes->ops->view            = SNESView_NCG;
449fef7b6d8SPeter Brune   snes->ops->solve           = SNESSolve_NCG;
450fef7b6d8SPeter Brune   snes->ops->reset           = SNESReset_NCG;
451fef7b6d8SPeter Brune 
452fef7b6d8SPeter Brune   snes->usesksp              = PETSC_FALSE;
453fef7b6d8SPeter Brune   snes->usespc               = PETSC_TRUE;
454fef7b6d8SPeter Brune 
4550e444f03SPeter Brune   snes->max_funcs = 30000;
4560e444f03SPeter Brune   snes->max_its   = 10000;
4570e444f03SPeter Brune 
458fef7b6d8SPeter Brune   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
459fef7b6d8SPeter Brune   snes->data = (void*) neP;
460dfb256c7SPeter Brune   neP->monitor = PETSC_NULL;
461*1867fe5bSPeter Brune   neP->betatype = 1;
46270d3b23bSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NCG",SNESLineSearchSetType_NCG);CHKERRQ(ierr);
463d2140ca3SPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
464fef7b6d8SPeter Brune   PetscFunctionReturn(0);
465fef7b6d8SPeter Brune }
466fef7b6d8SPeter Brune EXTERN_C_END
467