xref: /petsc/src/snes/impls/ncg/snesncg.c (revision fef7b6d8e45eee5b2327135157d7cc4a059fd396)
1*fef7b6d8SPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h>
2*fef7b6d8SPeter Brune 
3*fef7b6d8SPeter Brune #undef __FUNCT__
4*fef7b6d8SPeter Brune #define __FUNCT__ "SNESReset_NCG"
5*fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes)
6*fef7b6d8SPeter Brune {
7*fef7b6d8SPeter Brune   PetscErrorCode ierr;
8*fef7b6d8SPeter Brune 
9*fef7b6d8SPeter Brune   PetscFunctionBegin;
10*fef7b6d8SPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
11*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
12*fef7b6d8SPeter Brune }
13*fef7b6d8SPeter Brune 
14*fef7b6d8SPeter Brune /*
15*fef7b6d8SPeter Brune   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().
16*fef7b6d8SPeter Brune 
17*fef7b6d8SPeter Brune   Input Parameter:
18*fef7b6d8SPeter Brune . snes - the SNES context
19*fef7b6d8SPeter Brune 
20*fef7b6d8SPeter Brune   Application Interface Routine: SNESDestroy()
21*fef7b6d8SPeter Brune */
22*fef7b6d8SPeter Brune #undef __FUNCT__
23*fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG"
24*fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes)
25*fef7b6d8SPeter Brune {
26*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
27*fef7b6d8SPeter Brune   SNES_NCG *neP = (SNES_NCG*)snes->data;
28*fef7b6d8SPeter Brune 
29*fef7b6d8SPeter Brune   PetscFunctionBegin;
30*fef7b6d8SPeter Brune   ierr = SNESReset_NCG(snes);CHKERRQ(ierr);
31*fef7b6d8SPeter Brune   ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr);
32*fef7b6d8SPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
33*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
34*fef7b6d8SPeter Brune }
35*fef7b6d8SPeter Brune 
36*fef7b6d8SPeter Brune /*
37*fef7b6d8SPeter Brune    SNESSetUp_NCG - Sets up the internal data structures for the later use
38*fef7b6d8SPeter Brune    of the SNESNCG nonlinear solver.
39*fef7b6d8SPeter Brune 
40*fef7b6d8SPeter Brune    Input Parameters:
41*fef7b6d8SPeter Brune +  snes - the SNES context
42*fef7b6d8SPeter Brune -  x - the solution vector
43*fef7b6d8SPeter Brune 
44*fef7b6d8SPeter Brune    Application Interface Routine: SNESSetUp()
45*fef7b6d8SPeter Brune  */
46*fef7b6d8SPeter Brune #undef __FUNCT__
47*fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG"
48*fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes)
49*fef7b6d8SPeter Brune {
50*fef7b6d8SPeter Brune   PetscErrorCode ierr;
51*fef7b6d8SPeter Brune 
52*fef7b6d8SPeter Brune   PetscFunctionBegin;
53*fef7b6d8SPeter Brune   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
54*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
55*fef7b6d8SPeter Brune }
56*fef7b6d8SPeter Brune 
57*fef7b6d8SPeter Brune EXTERN_C_BEGIN
58*fef7b6d8SPeter Brune #undef __FUNCT__
59*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSetMonitor_NCG"
60*fef7b6d8SPeter Brune PetscErrorCode  SNESLineSearchSetMonitor_NCG(SNES snes,PetscBool  flg)
61*fef7b6d8SPeter Brune {
62*fef7b6d8SPeter Brune   SNES_NCG *neP = (SNES_NCG*)snes->data;
63*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
64*fef7b6d8SPeter Brune 
65*fef7b6d8SPeter Brune   PetscFunctionBegin;
66*fef7b6d8SPeter Brune   if (flg && !neP->monitor) {
67*fef7b6d8SPeter Brune     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&neP->monitor);CHKERRQ(ierr);
68*fef7b6d8SPeter Brune   } else if (!flg && neP->monitor) {
69*fef7b6d8SPeter Brune     ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr);
70*fef7b6d8SPeter Brune   }
71*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
72*fef7b6d8SPeter Brune }
73*fef7b6d8SPeter Brune EXTERN_C_END
74*fef7b6d8SPeter Brune 
75*fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNoNorms_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
76*fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNo_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
77*fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchQuadratic_NCG(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
78*fef7b6d8SPeter Brune /*
79*fef7b6d8SPeter Brune   SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method.
80*fef7b6d8SPeter Brune 
81*fef7b6d8SPeter Brune   Input Parameter:
82*fef7b6d8SPeter Brune . snes - the SNES context
83*fef7b6d8SPeter Brune 
84*fef7b6d8SPeter Brune   Application Interface Routine: SNESSetFromOptions()
85*fef7b6d8SPeter Brune */
86*fef7b6d8SPeter Brune #undef __FUNCT__
87*fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG"
88*fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
89*fef7b6d8SPeter Brune {
90*fef7b6d8SPeter Brune   SNES_NCG   *ls = (SNES_NCG *)snes->data;
91*fef7b6d8SPeter Brune   SNESLineSearchType indx;
92*fef7b6d8SPeter Brune   PetscBool          flg,set;
93*fef7b6d8SPeter Brune   PetscErrorCode     ierr;
94*fef7b6d8SPeter Brune 
95*fef7b6d8SPeter Brune   PetscFunctionBegin;
96*fef7b6d8SPeter Brune   ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr);
97*fef7b6d8SPeter Brune     ierr = PetscOptionsEnum("-snes_ls","NCG line search type","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)ls->type,(PetscEnum*)&indx,&flg);CHKERRQ(ierr);
98*fef7b6d8SPeter Brune     if (flg) {
99*fef7b6d8SPeter Brune       ls->type = indx;
100*fef7b6d8SPeter Brune       switch (indx) {
101*fef7b6d8SPeter Brune       case SNES_LS_BASIC:
102*fef7b6d8SPeter Brune         ierr = SNESLineSearchSet(snes,SNESLineSearchNo_NCG,PETSC_NULL);CHKERRQ(ierr);
103*fef7b6d8SPeter Brune         break;
104*fef7b6d8SPeter Brune       case SNES_LS_BASIC_NONORMS:
105*fef7b6d8SPeter Brune         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_NCG,PETSC_NULL);CHKERRQ(ierr);
106*fef7b6d8SPeter Brune         break;
107*fef7b6d8SPeter Brune       case SNES_LS_QUADRATIC:
108*fef7b6d8SPeter Brune         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NCG,PETSC_NULL);CHKERRQ(ierr);
109*fef7b6d8SPeter Brune         break;
110*fef7b6d8SPeter Brune       case SNES_LS_CUBIC:
111*fef7b6d8SPeter Brune         SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for cubic search");
112*fef7b6d8SPeter Brune         break;
113*fef7b6d8SPeter Brune       }
114*fef7b6d8SPeter Brune     }
115*fef7b6d8SPeter Brune     ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",ls->damping,&ls->damping,&flg);CHKERRQ(ierr);
116*fef7b6d8SPeter Brune     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
117*fef7b6d8SPeter Brune     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
118*fef7b6d8SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
119*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
120*fef7b6d8SPeter Brune }
121*fef7b6d8SPeter Brune 
122*fef7b6d8SPeter Brune /*
123*fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
124*fef7b6d8SPeter Brune 
125*fef7b6d8SPeter Brune   Input Parameters:
126*fef7b6d8SPeter Brune + SNES - the SNES context
127*fef7b6d8SPeter Brune - viewer - visualization context
128*fef7b6d8SPeter Brune 
129*fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
130*fef7b6d8SPeter Brune */
131*fef7b6d8SPeter Brune #undef __FUNCT__
132*fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
133*fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
134*fef7b6d8SPeter Brune {
135*fef7b6d8SPeter Brune   SNES_NCG *ls = (SNES_NCG *)snes->data;
136*fef7b6d8SPeter Brune   PetscBool        iascii;
137*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
138*fef7b6d8SPeter Brune 
139*fef7b6d8SPeter Brune   PetscFunctionBegin;
140*fef7b6d8SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
141*fef7b6d8SPeter Brune   if (iascii) {
142*fef7b6d8SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  line search type variant: %s\n", SNESLineSearchTypeName(ls->type));CHKERRQ(ierr);
143*fef7b6d8SPeter Brune   }
144*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
145*fef7b6d8SPeter Brune }
146*fef7b6d8SPeter Brune 
147*fef7b6d8SPeter Brune #undef __FUNCT__
148*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchNo_NCG"
149*fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNo_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
150*fef7b6d8SPeter Brune {
151*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
152*fef7b6d8SPeter Brune   SNES_NCG *neP = (SNES_NCG *) snes->data;
153*fef7b6d8SPeter Brune 
154*fef7b6d8SPeter Brune   PetscFunctionBegin;
155*fef7b6d8SPeter Brune   ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr);
156*fef7b6d8SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
157*fef7b6d8SPeter Brune   ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
158*fef7b6d8SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm");
159*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
160*fef7b6d8SPeter Brune }
161*fef7b6d8SPeter Brune 
162*fef7b6d8SPeter Brune #undef __FUNCT__
163*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchNoNorms_NCG"
164*fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchNoNorms_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
165*fef7b6d8SPeter Brune {
166*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
167*fef7b6d8SPeter Brune   SNES_NCG *neP = (SNES_NCG *) snes->data;
168*fef7b6d8SPeter Brune 
169*fef7b6d8SPeter Brune   PetscFunctionBegin;
170*fef7b6d8SPeter Brune   ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr);
171*fef7b6d8SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
172*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
173*fef7b6d8SPeter Brune }
174*fef7b6d8SPeter Brune 
175*fef7b6d8SPeter Brune #undef __FUNCT__
176*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_NCG"
177*fef7b6d8SPeter 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)
178*fef7b6d8SPeter Brune {
179*fef7b6d8SPeter Brune   PetscInt         i;
180*fef7b6d8SPeter Brune   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
181*fef7b6d8SPeter Brune   PetscReal        norms[3];
182*fef7b6d8SPeter Brune   PetscReal        alpha,a,b;
183*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
184*fef7b6d8SPeter Brune   SNES_NCG *neP = (SNES_NCG *) snes->data;
185*fef7b6d8SPeter Brune 
186*fef7b6d8SPeter Brune   PetscFunctionBegin;
187*fef7b6d8SPeter Brune   norms[0]  = fnorm;
188*fef7b6d8SPeter Brune   for(i=1; i < 3; ++i) {
189*fef7b6d8SPeter Brune     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
190*fef7b6d8SPeter Brune     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
191*fef7b6d8SPeter Brune     ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr);
192*fef7b6d8SPeter Brune   }
193*fef7b6d8SPeter Brune   for(i = 0; i < 3; ++i) {
194*fef7b6d8SPeter Brune     norms[i] = PetscSqr(norms[i]);
195*fef7b6d8SPeter Brune   }
196*fef7b6d8SPeter Brune   /* Fit a quadratic:
197*fef7b6d8SPeter Brune        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
198*fef7b6d8SPeter 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)
199*fef7b6d8SPeter 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)
200*fef7b6d8SPeter Brune        c = y_0
201*fef7b6d8SPeter Brune        x_min = -b/2a
202*fef7b6d8SPeter Brune 
203*fef7b6d8SPeter Brune        If we let x_{0,1,2} = 0, 0.5, 1.0
204*fef7b6d8SPeter Brune        a = 2 y_2 - 4 y_1 + 2 y_0
205*fef7b6d8SPeter Brune        b =  -y_2 + 4 y_1 - 3 y_0
206*fef7b6d8SPeter Brune        c =   y_0
207*fef7b6d8SPeter Brune   */
208*fef7b6d8SPeter 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]));
209*fef7b6d8SPeter 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]);
210*fef7b6d8SPeter Brune   /* Check for positive a (concave up) */
211*fef7b6d8SPeter Brune   if (a >= 0.0) {
212*fef7b6d8SPeter Brune     alpha = -b/(2.0*a);
213*fef7b6d8SPeter Brune     alpha = PetscMin(alpha, alphas[2]);
214*fef7b6d8SPeter Brune     alpha = PetscMax(alpha, alphas[0]);
215*fef7b6d8SPeter Brune   } else {
216*fef7b6d8SPeter Brune     alpha = 1.0;
217*fef7b6d8SPeter Brune   }
218*fef7b6d8SPeter Brune   if (neP->monitor) {
219*fef7b6d8SPeter Brune     ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
220*fef7b6d8SPeter Brune     ierr = PetscViewerASCIIPrintf(neP->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);
221*fef7b6d8SPeter Brune     ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
222*fef7b6d8SPeter Brune   }
223*fef7b6d8SPeter Brune   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
224*fef7b6d8SPeter Brune   if (alpha != 1.0) {
225*fef7b6d8SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
226*fef7b6d8SPeter Brune     ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
227*fef7b6d8SPeter Brune   } else {
228*fef7b6d8SPeter Brune     *gnorm = PetscSqrtReal(norms[2]);
229*fef7b6d8SPeter Brune   }
230*fef7b6d8SPeter Brune   if (alpha == 0.0) *flag = PETSC_FALSE;
231*fef7b6d8SPeter Brune   else              *flag = PETSC_TRUE;
232*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
233*fef7b6d8SPeter Brune }
234*fef7b6d8SPeter Brune 
235*fef7b6d8SPeter Brune /*
236*fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
237*fef7b6d8SPeter Brune 
238*fef7b6d8SPeter Brune   Input Parameters:
239*fef7b6d8SPeter Brune . snes - the SNES context
240*fef7b6d8SPeter Brune 
241*fef7b6d8SPeter Brune   Output Parameter:
242*fef7b6d8SPeter Brune . outits - number of iterations until termination
243*fef7b6d8SPeter Brune 
244*fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
245*fef7b6d8SPeter Brune */
246*fef7b6d8SPeter Brune #undef __FUNCT__
247*fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
248*fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
249*fef7b6d8SPeter Brune {
250*fef7b6d8SPeter Brune   SNES_NCG            *neP = (SNES_NCG *) snes->data;
251*fef7b6d8SPeter Brune   Vec                 X, dX, lX, F, W;
252*fef7b6d8SPeter Brune   PetscReal           fnorm, dummyNorm;
253*fef7b6d8SPeter Brune   PetscScalar         dXdot, dXdot_old;
254*fef7b6d8SPeter Brune   PetscInt            maxits, i;
255*fef7b6d8SPeter Brune   PetscErrorCode      ierr;
256*fef7b6d8SPeter Brune   SNESConvergedReason reason;
257*fef7b6d8SPeter Brune   PetscBool           lsSuccess = PETSC_TRUE;
258*fef7b6d8SPeter Brune   PetscFunctionBegin;
259*fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
260*fef7b6d8SPeter Brune 
261*fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
262*fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
263*fef7b6d8SPeter Brune   dX     = snes->vec_sol_update;     /* the steepest direction */
264*fef7b6d8SPeter Brune   lX     = snes->work[1];            /* the conjugate direction */
265*fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
266*fef7b6d8SPeter Brune   W      = snes->work[0];            /* work vector for the line search */
267*fef7b6d8SPeter Brune 
268*fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
269*fef7b6d8SPeter Brune   snes->iter = 0;
270*fef7b6d8SPeter Brune   snes->norm = 0.;
271*fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
272*fef7b6d8SPeter Brune 
273*fef7b6d8SPeter Brune   /* compute the initial function and preconditioned update delX */
274*fef7b6d8SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
275*fef7b6d8SPeter Brune   if (snes->domainerror) {
276*fef7b6d8SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
277*fef7b6d8SPeter Brune     PetscFunctionReturn(0);
278*fef7b6d8SPeter Brune   }
279*fef7b6d8SPeter Brune   /* convergence test */
280*fef7b6d8SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
281*fef7b6d8SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
282*fef7b6d8SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
283*fef7b6d8SPeter Brune   snes->norm = fnorm;
284*fef7b6d8SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
285*fef7b6d8SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
286*fef7b6d8SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
287*fef7b6d8SPeter Brune 
288*fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
289*fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
290*fef7b6d8SPeter Brune   /* test convergence */
291*fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
292*fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
293*fef7b6d8SPeter Brune 
294*fef7b6d8SPeter Brune   /* Call general purpose update function */
295*fef7b6d8SPeter Brune   if (snes->ops->update) {
296*fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
297*fef7b6d8SPeter Brune   }
298*fef7b6d8SPeter Brune 
299*fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
300*fef7b6d8SPeter Brune   if (!snes->pc) {
301*fef7b6d8SPeter Brune     ierr = VecCopy(F, lX);CHKERRQ(ierr);
302*fef7b6d8SPeter Brune     ierr = VecScale(lX, -1.0);CHKERRQ(ierr);
303*fef7b6d8SPeter Brune   } else {
304*fef7b6d8SPeter Brune     ierr = VecCopy(X, lX);CHKERRQ(ierr);
305*fef7b6d8SPeter Brune     ierr = SNESSolve(snes->pc, 0, lX);CHKERRQ(ierr);
306*fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
307*fef7b6d8SPeter Brune     ierr = VecAXPY(dX, -1.0, X);CHKERRQ(ierr);
308*fef7b6d8SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
309*fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
310*fef7b6d8SPeter Brune       PetscFunctionReturn(0);
311*fef7b6d8SPeter Brune     }
312*fef7b6d8SPeter Brune   }
313*fef7b6d8SPeter Brune   ierr = VecDot(lX, lX, &dXdot);CHKERRQ(ierr);
314*fef7b6d8SPeter Brune   /* line search for the update */
315*fef7b6d8SPeter Brune   ierr = (*neP->LineSearch)(snes, neP->lsP, X, F, dX, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
316*fef7b6d8SPeter Brune 
317*fef7b6d8SPeter Brune   for(i = 1; i < maxits; i++) {
318*fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
319*fef7b6d8SPeter Brune     /* Call general purpose update function */
320*fef7b6d8SPeter Brune     if (snes->ops->update) {
321*fef7b6d8SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
322*fef7b6d8SPeter Brune     }
323*fef7b6d8SPeter Brune     if (!snes->pc) {
324*fef7b6d8SPeter Brune       ierr = VecCopy(F,dX);CHKERRQ(ierr);
325*fef7b6d8SPeter Brune       ierr = VecScale(dX,-1.0);CHKERRQ(ierr);
326*fef7b6d8SPeter Brune     } else {
327*fef7b6d8SPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
328*fef7b6d8SPeter Brune       ierr = SNESSolve(snes->pc, 0, dX);CHKERRQ(ierr);
329*fef7b6d8SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
330*fef7b6d8SPeter Brune       if (reason < 0) {
331*fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
332*fef7b6d8SPeter Brune         PetscFunctionReturn(0);
333*fef7b6d8SPeter Brune       }
334*fef7b6d8SPeter Brune       ierr = VecAXPY(dX,-1.0,X);CHKERRQ(ierr);
335*fef7b6d8SPeter Brune     }
336*fef7b6d8SPeter Brune 
337*fef7b6d8SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = (dX, dX) / (DX_old, dX_old) (Fletcher-Reeves Update)*/
338*fef7b6d8SPeter Brune     dXdot_old = dXdot;
339*fef7b6d8SPeter Brune     ierr = VecDot(dX, dX, &dXdot);CHKERRQ(ierr);
340*fef7b6d8SPeter Brune     ierr = VecAXPY(lX, dXdot / dXdot_old, dX);CHKERRQ(ierr);
341*fef7b6d8SPeter Brune     /* line search for the proper contribution of lX to the solution */
342*fef7b6d8SPeter Brune     ierr = (*neP->LineSearch)(snes, neP->lsP, X, F, lX, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
343*fef7b6d8SPeter Brune     if (!lsSuccess) {
344*fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
345*fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
346*fef7b6d8SPeter Brune         break;
347*fef7b6d8SPeter Brune       }
348*fef7b6d8SPeter Brune     }
349*fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
350*fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
351*fef7b6d8SPeter Brune       break;
352*fef7b6d8SPeter Brune     }
353*fef7b6d8SPeter Brune     if (snes->domainerror) {
354*fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
355*fef7b6d8SPeter Brune       PetscFunctionReturn(0);
356*fef7b6d8SPeter Brune     }
357*fef7b6d8SPeter Brune     /* Monitor convergence */
358*fef7b6d8SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
359*fef7b6d8SPeter Brune     snes->iter = i+1;
360*fef7b6d8SPeter Brune     snes->norm = fnorm;
361*fef7b6d8SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
362*fef7b6d8SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
363*fef7b6d8SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
364*fef7b6d8SPeter Brune     /* Test for convergence */
365*fef7b6d8SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
366*fef7b6d8SPeter Brune     if (snes->reason) break;
367*fef7b6d8SPeter Brune   }
368*fef7b6d8SPeter Brune   if (i == maxits) {
369*fef7b6d8SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
370*fef7b6d8SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
371*fef7b6d8SPeter Brune   }
372*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
373*fef7b6d8SPeter Brune }
374*fef7b6d8SPeter Brune 
375*fef7b6d8SPeter Brune typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,void*,PetscBool *);                 /* force argument to next function to not be extern C*/
376*fef7b6d8SPeter Brune EXTERN_C_BEGIN
377*fef7b6d8SPeter Brune #undef __FUNCT__
378*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSetPreCheck_NCG"
379*fef7b6d8SPeter Brune PetscErrorCode  SNESLineSearchSetPreCheck_NCG(SNES snes, FCN1 func, void *checkctx)
380*fef7b6d8SPeter Brune {
381*fef7b6d8SPeter Brune   PetscFunctionBegin;
382*fef7b6d8SPeter Brune   ((SNES_NCG *)(snes->data))->precheckstep = func;
383*fef7b6d8SPeter Brune   ((SNES_NCG *)(snes->data))->precheck     = checkctx;
384*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
385*fef7b6d8SPeter Brune }
386*fef7b6d8SPeter Brune EXTERN_C_END
387*fef7b6d8SPeter Brune 
388*fef7b6d8SPeter Brune typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
389*fef7b6d8SPeter Brune EXTERN_C_BEGIN
390*fef7b6d8SPeter Brune #undef __FUNCT__
391*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSet_NCG"
392*fef7b6d8SPeter Brune PetscErrorCode  SNESLineSearchSet_NCG(SNES snes, FCN2 func, void *lsctx)
393*fef7b6d8SPeter Brune {
394*fef7b6d8SPeter Brune   PetscFunctionBegin;
395*fef7b6d8SPeter Brune   ((SNES_NCG *)(snes->data))->LineSearch = func;
396*fef7b6d8SPeter Brune   ((SNES_NCG *)(snes->data))->lsP        = lsctx;
397*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
398*fef7b6d8SPeter Brune }
399*fef7b6d8SPeter Brune EXTERN_C_END
400*fef7b6d8SPeter Brune 
401*fef7b6d8SPeter Brune typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/
402*fef7b6d8SPeter Brune EXTERN_C_BEGIN
403*fef7b6d8SPeter Brune #undef __FUNCT__
404*fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchSetPostCheck_NCG"
405*fef7b6d8SPeter Brune PetscErrorCode  SNESLineSearchSetPostCheck_NCG(SNES snes, FCN3 func, void *checkctx)
406*fef7b6d8SPeter Brune {
407*fef7b6d8SPeter Brune   PetscFunctionBegin;
408*fef7b6d8SPeter Brune   ((SNES_NCG *)(snes->data))->postcheckstep = func;
409*fef7b6d8SPeter Brune   ((SNES_NCG *)(snes->data))->postcheck     = checkctx;
410*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
411*fef7b6d8SPeter Brune }
412*fef7b6d8SPeter Brune EXTERN_C_END
413*fef7b6d8SPeter Brune 
414*fef7b6d8SPeter Brune /*MC
415*fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
416*fef7b6d8SPeter Brune 
417*fef7b6d8SPeter Brune   Level: beginner
418*fef7b6d8SPeter Brune 
419*fef7b6d8SPeter Brune   Options Database:
420*fef7b6d8SPeter Brune +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
421*fef7b6d8SPeter Brune -   -snes_ls <basic,basicnormnorms,quadratic>
422*fef7b6d8SPeter Brune 
423*fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
424*fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
425*fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
426*fef7b6d8SPeter Brune 
427*fef7b6d8SPeter Brune 
428*fef7b6d8SPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
429*fef7b6d8SPeter Brune M*/
430*fef7b6d8SPeter Brune EXTERN_C_BEGIN
431*fef7b6d8SPeter Brune #undef __FUNCT__
432*fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
433*fef7b6d8SPeter Brune PetscErrorCode  SNESCreate_NCG(SNES snes)
434*fef7b6d8SPeter Brune {
435*fef7b6d8SPeter Brune   SNES_NCG *neP;
436*fef7b6d8SPeter Brune   PetscErrorCode   ierr;
437*fef7b6d8SPeter Brune 
438*fef7b6d8SPeter Brune   PetscFunctionBegin;
439*fef7b6d8SPeter Brune   snes->ops->destroy         = SNESDestroy_NCG;
440*fef7b6d8SPeter Brune   snes->ops->setup           = SNESSetUp_NCG;
441*fef7b6d8SPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_NCG;
442*fef7b6d8SPeter Brune   snes->ops->view            = SNESView_NCG;
443*fef7b6d8SPeter Brune   snes->ops->solve           = SNESSolve_NCG;
444*fef7b6d8SPeter Brune   snes->ops->reset           = SNESReset_NCG;
445*fef7b6d8SPeter Brune 
446*fef7b6d8SPeter Brune   snes->usesksp              = PETSC_FALSE;
447*fef7b6d8SPeter Brune   snes->usespc               = PETSC_TRUE;
448*fef7b6d8SPeter Brune 
449*fef7b6d8SPeter Brune   ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
450*fef7b6d8SPeter Brune   snes->data = (void*) neP;
451*fef7b6d8SPeter Brune   neP->maxstep       = 1.e8;
452*fef7b6d8SPeter Brune   neP->steptol       = 1.e-12;
453*fef7b6d8SPeter Brune   neP->type          = SNES_LS_QUADRATIC;
454*fef7b6d8SPeter Brune   neP->damping       = 1.0;
455*fef7b6d8SPeter Brune   neP->LineSearch    = SNESLineSearchQuadratic_NCG;
456*fef7b6d8SPeter Brune   neP->lsP           = PETSC_NULL;
457*fef7b6d8SPeter Brune   neP->postcheckstep = PETSC_NULL;
458*fef7b6d8SPeter Brune   neP->postcheck     = PETSC_NULL;
459*fef7b6d8SPeter Brune   neP->precheckstep  = PETSC_NULL;
460*fef7b6d8SPeter Brune   neP->precheck      = PETSC_NULL;
461*fef7b6d8SPeter Brune 
462*fef7b6d8SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_NCG",SNESLineSearchSetMonitor_NCG);CHKERRQ(ierr);
463*fef7b6d8SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_NCG",SNESLineSearchSet_NCG);CHKERRQ(ierr);
464*fef7b6d8SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_NCG",SNESLineSearchSetPostCheck_NCG);CHKERRQ(ierr);
465*fef7b6d8SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_NCG",SNESLineSearchSetPreCheck_NCG);CHKERRQ(ierr);
466*fef7b6d8SPeter Brune   PetscFunctionReturn(0);
467*fef7b6d8SPeter Brune }
468*fef7b6d8SPeter Brune EXTERN_C_END
469