xref: /petsc/src/snes/impls/ls/ls.c (revision a5892487d3a72298dc86de4f4aaa34ef49eae3db)
1 
2 #include <../src/snes/impls/ls/lsimpl.h>
3 
4 /*
5      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
6     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
7     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
8     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
9 */
10 #undef __FUNCT__
11 #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12 PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
13 {
14   PetscReal      a1;
15   PetscErrorCode ierr;
16   PetscBool      hastranspose;
17 
18   PetscFunctionBegin;
19   *ismin = PETSC_FALSE;
20   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
21   if (hastranspose) {
22     /* Compute || J^T F|| */
23     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
24     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
25     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
26     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
27   } else {
28     Vec         work;
29     PetscScalar result;
30     PetscReal   wnorm;
31 
32     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
33     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
34     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
35     ierr = MatMult(A,W,work);CHKERRQ(ierr);
36     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
37     ierr = VecDestroy(&work);CHKERRQ(ierr);
38     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
39     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
40     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46      Checks if J^T(F - J*X) = 0
47 */
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESLSCheckResidual_Private"
50 PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
51 {
52   PetscReal      a1,a2;
53   PetscErrorCode ierr;
54   PetscBool      hastranspose;
55 
56   PetscFunctionBegin;
57   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
58   if (hastranspose) {
59     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
60     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
61 
62     /* Compute || J^T W|| */
63     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
64     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
65     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
66     if (a1 != 0.0) {
67       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
68     }
69   }
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "SNESLineSearchCubic_LS"
75 /*@C
76    SNESLineSearchCubic - Performs a cubic line search (default line search method).
77 
78    Collective on SNES
79 
80    Input Parameters:
81 +  snes - nonlinear context
82 .  lsctx - optional context for line search (not used here)
83 .  x - current iterate
84 .  f - residual evaluated at x
85 .  y - search direction
86 .  fnorm - 2-norm of f
87 -  xnorm - norm of x if known, otherwise 0
88 
89    Output Parameters:
90 +  g - residual evaluated at new iterate y
91 .  w - new iterate
92 .  gnorm - 2-norm of g
93 .  ynorm - 2-norm of search length
94 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
95 
96    Options Database Key:
97 +  -snes_ls cubic - Activates SNESLineSearchCubic()
98 .   -snes_ls_alpha <alpha> - Sets alpha
99 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
100 -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
101 
102 
103    Notes:
104    This line search is taken from "Numerical Methods for Unconstrained
105    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
106 
107    Level: advanced
108 
109 .keywords: SNES, nonlinear, line search, cubic
110 
111 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
112 @*/
113 PetscErrorCode  SNESLineSearchCubic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
114 {
115   /*
116      Note that for line search purposes we work with with the related
117      minimization problem:
118         min  z(x):  R^n -> R,
119      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
120    */
121 
122   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
123   PetscReal      minlambda,lambda,lambdatemp;
124 #if defined(PETSC_USE_COMPLEX)
125   PetscScalar    cinitslope;
126 #endif
127   PetscErrorCode ierr;
128   PetscInt       count;
129   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
130   MPI_Comm       comm;
131 
132   PetscFunctionBegin;
133   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
134   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
135   *flag   = PETSC_TRUE;
136 
137   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
138   if (*ynorm == 0.0) {
139     if (snes->ls_monitor) {
140       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
141       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
142       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
143     }
144     *gnorm = fnorm;
145     ierr   = VecCopy(x,w);CHKERRQ(ierr);
146     ierr   = VecCopy(f,g);CHKERRQ(ierr);
147     *flag  = PETSC_FALSE;
148     goto theend1;
149   }
150   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
151     if (snes->ls_monitor) {
152       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
153       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
154       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
155     }
156     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
157     *ynorm = snes->maxstep;
158   }
159   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
160   minlambda = snes->steptol/rellength;
161   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
162 #if defined(PETSC_USE_COMPLEX)
163   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
164   initslope = PetscRealPart(cinitslope);
165 #else
166   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
167 #endif
168   if (initslope > 0.0)  initslope = -initslope;
169   if (initslope == 0.0) initslope = -1.0;
170 
171   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
172   if (snes->nfuncs >= snes->max_funcs) {
173     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
174     *flag = PETSC_FALSE;
175     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
176     goto theend1;
177   }
178   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
179   if (snes->domainerror) {
180     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
181     PetscFunctionReturn(0);
182   }
183   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
184   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
185   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
186   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
187     if (snes->ls_monitor) {
188       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
189       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
190       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
191     }
192     goto theend1;
193   }
194 
195   /* Fit points with quadratic */
196   lambda     = 1.0;
197   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
198   lambdaprev = lambda;
199   gnormprev  = *gnorm;
200   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
201   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
202   else                         lambda = lambdatemp;
203 
204   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
205   if (snes->nfuncs >= snes->max_funcs) {
206     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
207     *flag = PETSC_FALSE;
208     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
209     goto theend1;
210   }
211   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
212   if (snes->domainerror) {
213     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
214     PetscFunctionReturn(0);
215   }
216   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
217   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
218   if (snes->ls_monitor) {
219     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
220     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr);
221     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
222   }
223   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
224     if (snes->ls_monitor) {
225       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
226       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
227       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
228     }
229     goto theend1;
230   }
231 
232   /* Fit points with cubic */
233   count = 1;
234   while (PETSC_TRUE) {
235     if (lambda <= minlambda) {
236       if (snes->ls_monitor) {
237         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
238 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
239 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr);
240         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
241       }
242       *flag = PETSC_FALSE;
243       break;
244     }
245     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
246     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
247     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
248     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
249     d  = b*b - 3*a*initslope;
250     if (d < 0.0) d = 0.0;
251     if (a == 0.0) {
252       lambdatemp = -initslope/(2.0*b);
253     } else {
254       lambdatemp = (-b + sqrt(d))/(3.0*a);
255     }
256     lambdaprev = lambda;
257     gnormprev  = *gnorm;
258     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
259     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
260     else                         lambda     = lambdatemp;
261     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
262     if (snes->nfuncs >= snes->max_funcs) {
263       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
264       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
265       *flag = PETSC_FALSE;
266       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
267       break;
268     }
269     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
270     if (snes->domainerror) {
271       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
272       PetscFunctionReturn(0);
273     }
274     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
275     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
276     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* is reduction enough? */
277       if (snes->ls_monitor) {
278 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
279       }
280       break;
281     } else {
282       if (snes->ls_monitor) {
283         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
284       }
285     }
286     count++;
287   }
288   theend1:
289   /* Optional user-defined check for line search step validity */
290   if (snes->ops->postcheckstep && *flag) {
291     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
292     if (changed_y) {
293       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
294     }
295     if (changed_y || changed_w) { /* recompute the function if the step has changed */
296       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
297       if (snes->domainerror) {
298         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
299         PetscFunctionReturn(0);
300       }
301       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
302       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
303       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
304       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
305       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
306     }
307   }
308   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
309   PetscFunctionReturn(0);
310 }
311 /* -------------------------------------------------------------------------- */
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "SNESLineSearchQuadratic_LS"
315 /*@C
316    SNESLineSearchQuadratic_LS - Performs a quadratic line search.
317 
318    Collective on SNES and Vec
319 
320    Input Parameters:
321 +  snes - the SNES context
322 .  lsctx - optional context for line search (not used here)
323 .  x - current iterate
324 .  f - residual evaluated at x
325 .  y - search direction
326 .  fnorm - 2-norm of f
327 -  xnorm - norm of x if known, otherwise 0
328 
329    Output Parameters:
330 +  g - residual evaluated at new iterate w
331 .  w - new iterate (x + lambda*y)
332 .  gnorm - 2-norm of g
333 .  ynorm - 2-norm of search length
334 -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
335 
336    Options Database Keys:
337 +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
338 .   -snes_ls_alpha <alpha> - Sets alpha
339 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
340 -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
341 
342    Notes:
343    Use SNESLineSearchSet() to set this routine within the SNESLS method.
344 
345    Level: advanced
346 
347 .keywords: SNES, nonlinear, quadratic, line search
348 
349 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
350 @*/
351 PetscErrorCode  SNESLineSearchQuadratic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
352 {
353   /*
354      Note that for line search purposes we work with with the related
355      minimization problem:
356         min  z(x):  R^n -> R,
357      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
358    */
359   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
360 #if defined(PETSC_USE_COMPLEX)
361   PetscScalar    cinitslope;
362 #endif
363   PetscErrorCode ierr;
364   PetscInt       count;
365   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
366 
367   PetscFunctionBegin;
368   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
369   *flag   = PETSC_TRUE;
370 
371   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
372   if (*ynorm == 0.0) {
373     if (snes->ls_monitor) {
374       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
375       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
376       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
377     }
378     *gnorm = fnorm;
379     ierr   = VecCopy(x,w);CHKERRQ(ierr);
380     ierr   = VecCopy(f,g);CHKERRQ(ierr);
381     *flag  = PETSC_FALSE;
382     goto theend2;
383   }
384   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
385     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
386     *ynorm = snes->maxstep;
387   }
388   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
389   minlambda = snes->steptol/rellength;
390   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
391 #if defined(PETSC_USE_COMPLEX)
392   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
393   initslope = PetscRealPart(cinitslope);
394 #else
395   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
396 #endif
397   if (initslope > 0.0)  initslope = -initslope;
398   if (initslope == 0.0) initslope = -1.0;
399   ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr);
400 
401   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
402   if (snes->nfuncs >= snes->max_funcs) {
403     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
404     *flag = PETSC_FALSE;
405     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
406     goto theend2;
407   }
408   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
409   if (snes->domainerror) {
410     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
411     PetscFunctionReturn(0);
412   }
413   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
414   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
415   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
416     if (snes->ls_monitor) {
417       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
418       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
419       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
420     }
421     goto theend2;
422   }
423 
424   /* Fit points with quadratic */
425   lambda = 1.0;
426   count = 1;
427   while (PETSC_TRUE) {
428     if (lambda <= minlambda) { /* bad luck; use full step */
429       if (snes->ls_monitor) {
430         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
431         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
432         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%14.12e, gnorm=%14.12e, ynorm=%14.12e, lambda=%14.12e, initial slope=%14.12e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
433         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
434       }
435       ierr = VecCopy(x,w);CHKERRQ(ierr);
436       *flag = PETSC_FALSE;
437       break;
438     }
439     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
440     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
441     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
442     else                         lambda     = lambdatemp;
443 
444     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
445     if (snes->nfuncs >= snes->max_funcs) {
446       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
447       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
448       *flag = PETSC_FALSE;
449       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
450       break;
451     }
452     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
453     if (snes->domainerror) {
454       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
455       PetscFunctionReturn(0);
456     }
457     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
458     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
459     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
460       if (snes->ls_monitor) {
461         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
462         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr);
463         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
464       }
465       break;
466     }
467     count++;
468   }
469   theend2:
470   /* Optional user-defined check for line search step validity */
471   if (snes->ops->postcheckstep) {
472     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
473     if (changed_y) {
474       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
475     }
476     if (changed_y || changed_w) { /* recompute the function if the step has changed */
477       ierr = SNESComputeFunction(snes,w,g);
478       if (snes->domainerror) {
479         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
480         PetscFunctionReturn(0);
481       }
482       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
483       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
484       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
485       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
486       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
487     }
488   }
489   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
490   PetscFunctionReturn(0);
491 }
492 
493 
494 /*  --------------------------------------------------------------------
495 
496      This file implements a truncated Newton method with a line search,
497      for solving a system of nonlinear equations, using the KSP, Vec,
498      and Mat interfaces for linear solvers, vectors, and matrices,
499      respectively.
500 
501      The following basic routines are required for each nonlinear solver:
502           SNESCreate_XXX()          - Creates a nonlinear solver context
503           SNESSetFromOptions_XXX()  - Sets runtime options
504           SNESSolve_XXX()           - Solves the nonlinear system
505           SNESDestroy_XXX()         - Destroys the nonlinear solver context
506      The suffix "_XXX" denotes a particular implementation, in this case
507      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
508      systems of nonlinear equations with a line search (LS) method.
509      These routines are actually called via the common user interface
510      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
511      SNESDestroy(), so the application code interface remains identical
512      for all nonlinear solvers.
513 
514      Another key routine is:
515           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
516      by setting data structures and options.   The interface routine SNESSetUp()
517      is not usually called directly by the user, but instead is called by
518      SNESSolve() if necessary.
519 
520      Additional basic routines are:
521           SNESView_XXX()            - Prints details of runtime options that
522                                       have actually been used.
523      These are called by application codes via the interface routines
524      SNESView().
525 
526      The various types of solvers (preconditioners, Krylov subspace methods,
527      nonlinear solvers, timesteppers) are all organized similarly, so the
528      above description applies to these categories also.
529 
530     -------------------------------------------------------------------- */
531 /*
532    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
533    method with a line search.
534 
535    Input Parameters:
536 .  snes - the SNES context
537 
538    Output Parameter:
539 .  outits - number of iterations until termination
540 
541    Application Interface Routine: SNESSolve()
542 
543    Notes:
544    This implements essentially a truncated Newton method with a
545    line search.  By default a cubic backtracking line search
546    is employed, as described in the text "Numerical Methods for
547    Unconstrained Optimization and Nonlinear Equations" by Dennis
548    and Schnabel.
549 */
550 #undef __FUNCT__
551 #define __FUNCT__ "SNESSolve_LS"
552 PetscErrorCode SNESSolve_LS(SNES snes)
553 {
554   PetscErrorCode     ierr;
555   PetscInt           maxits,i,lits;
556   PetscBool          lssucceed;
557   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
558   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
559   Vec                Y,X,F,G,W;
560   KSPConvergedReason kspreason;
561 
562   PetscFunctionBegin;
563   snes->numFailures            = 0;
564   snes->numLinearSolveFailures = 0;
565   snes->reason                 = SNES_CONVERGED_ITERATING;
566 
567   maxits	= snes->max_its;	/* maximum number of iterations */
568   X		= snes->vec_sol;	/* solution vector */
569   F		= snes->vec_func;	/* residual vector */
570   Y		= snes->work[0];	/* work vectors */
571   G		= snes->work[1];
572   W		= snes->work[2];
573 
574   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
575   snes->iter = 0;
576   snes->norm = 0.0;
577   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
578   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
579   if (snes->domainerror) {
580     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
581     PetscFunctionReturn(0);
582   }
583   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
584   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
585   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
586   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
587   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
588   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
589   snes->norm = fnorm;
590   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
591   SNESLogConvHistory(snes,fnorm,0);
592   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
593 
594   /* set parameter for default relative tolerance convergence test */
595   snes->ttol = fnorm*snes->rtol;
596   /* test convergence */
597   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
598   if (snes->reason) PetscFunctionReturn(0);
599 
600   for (i=0; i<maxits; i++) {
601 
602     /* Call general purpose update function */
603     if (snes->ops->update) {
604       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
605     }
606 
607     /* Solve J Y = F, where J is Jacobian matrix */
608     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
609     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
610     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
611     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
612     if (kspreason < 0) {
613       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
614         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
615         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
616         break;
617       }
618     }
619     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
620     snes->linear_its += lits;
621     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
622 
623     if (snes->ops->precheckstep) {
624       PetscBool  changed_y = PETSC_FALSE;
625       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
626     }
627 
628     if (PetscLogPrintInfo){
629       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
630     }
631 
632     /* Compute a (scaled) negative update in the line search routine:
633          Y <- X - lambda*Y
634        and evaluate G = function(Y) (depends on the line search).
635     */
636     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
637     ynorm = 1; gnorm = fnorm;
638     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
639     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
640     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
641     if (snes->domainerror) {
642       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
643       PetscFunctionReturn(0);
644     }
645     if (!lssucceed) {
646       if (++snes->numFailures >= snes->maxFailures) {
647 	PetscBool  ismin;
648         snes->reason = SNES_DIVERGED_LINE_SEARCH;
649         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
650         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
651         break;
652       }
653     }
654     /* Update function and solution vectors */
655     fnorm = gnorm;
656     ierr = VecCopy(G,F);CHKERRQ(ierr);
657     ierr = VecCopy(W,X);CHKERRQ(ierr);
658     /* Monitor convergence */
659     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
660     snes->iter = i+1;
661     snes->norm = fnorm;
662     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
663     SNESLogConvHistory(snes,snes->norm,lits);
664     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
665     /* Test for convergence, xnorm = || X || */
666     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
667     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
668     if (snes->reason) break;
669   }
670   if (i == maxits) {
671     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
672     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
673   }
674   PetscFunctionReturn(0);
675 }
676 /* -------------------------------------------------------------------------- */
677 /*
678    SNESSetUp_LS - Sets up the internal data structures for the later use
679    of the SNESLS nonlinear solver.
680 
681    Input Parameter:
682 .  snes - the SNES context
683 .  x - the solution vector
684 
685    Application Interface Routine: SNESSetUp()
686 
687    Notes:
688    For basic use of the SNES solvers, the user need not explicitly call
689    SNESSetUp(), since these actions will automatically occur during
690    the call to SNESSolve().
691  */
692 #undef __FUNCT__
693 #define __FUNCT__ "SNESSetUp_LS"
694 PetscErrorCode SNESSetUp_LS(SNES snes)
695 {
696   PetscErrorCode ierr;
697 
698   PetscFunctionBegin;
699   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
700   PetscFunctionReturn(0);
701 }
702 /* -------------------------------------------------------------------------- */
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "SNESReset_LS"
706 PetscErrorCode SNESReset_LS(SNES snes)
707 {
708   PetscErrorCode ierr;
709 
710   PetscFunctionBegin;
711   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
712   PetscFunctionReturn(0);
713 }
714 
715 /*
716    SNESDestroy_LS - Destroys the private SNES_LS context that was created
717    with SNESCreate_LS().
718 
719    Input Parameter:
720 .  snes - the SNES context
721 
722    Application Interface Routine: SNESDestroy()
723  */
724 #undef __FUNCT__
725 #define __FUNCT__ "SNESDestroy_LS"
726 PetscErrorCode SNESDestroy_LS(SNES snes)
727 {
728   PetscErrorCode ierr;
729 
730   PetscFunctionBegin;
731   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
732   ierr = PetscFree(snes->data);CHKERRQ(ierr);
733 
734   PetscFunctionReturn(0);
735 }
736 /* -------------------------------------------------------------------------- */
737 
738 /*
739    SNESView_LS - Prints info from the SNESLS data structure.
740 
741    Input Parameters:
742 .  SNES - the SNES context
743 .  viewer - visualization context
744 
745    Application Interface Routine: SNESView()
746 */
747 #undef __FUNCT__
748 #define __FUNCT__ "SNESView_LS"
749 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
750 {
751   const char     *cstr;
752   PetscErrorCode ierr;
753   PetscBool      iascii;
754 
755   PetscFunctionBegin;
756   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
757   if (iascii) {
758     cstr = SNESLineSearchTypeName(snes->ls_type);
759     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
760     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)snes->ls_alpha,(double)snes->maxstep,(double)snes->steptol);CHKERRQ(ierr);
761     ierr = PetscViewerASCIIPrintf(viewer,"  damping factor=%14.12e\n",(double)snes->damping);CHKERRQ(ierr);
762   }
763   PetscFunctionReturn(0);
764 }
765 /* -------------------------------------------------------------------------- */
766 /*
767    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
768 
769    Input Parameter:
770 .  snes - the SNES context
771 
772    Application Interface Routine: SNESSetFromOptions()
773 */
774 #undef __FUNCT__
775 #define __FUNCT__ "SNESSetFromOptions_LS"
776 static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
777 {
778   PetscErrorCode     ierr;
779 
780   PetscFunctionBegin;
781   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
782   ierr = PetscOptionsTail();CHKERRQ(ierr);
783   PetscFunctionReturn(0);
784 }
785 
786 EXTERN_C_BEGIN
787 extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
788 EXTERN_C_END
789 
790 /* -------------------------------------------------------------------------- */
791 /*MC
792       SNESLS - Newton based nonlinear solver that uses a line search
793 
794    Options Database:
795 +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
796 .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
797 .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
798 .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
799 .   -snes_ls_monitor - print information about progress of line searches
800 -   -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms
801 
802 
803     Notes: This is the default nonlinear solver in SNES
804 
805    Level: beginner
806 
807 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
808            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
809           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
810 
811 M*/
812 EXTERN_C_BEGIN
813 #undef __FUNCT__
814 #define __FUNCT__ "SNESCreate_LS"
815 PetscErrorCode  SNESCreate_LS(SNES snes)
816 {
817   PetscErrorCode ierr;
818   SNES_LS        *neP;
819 
820   PetscFunctionBegin;
821   snes->ops->setup           = SNESSetUp_LS;
822   snes->ops->solve           = SNESSolve_LS;
823   snes->ops->destroy         = SNESDestroy_LS;
824   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
825   snes->ops->view            = SNESView_LS;
826   snes->ops->reset           = SNESReset_LS;
827 
828   snes->usesksp                      = PETSC_TRUE;
829   snes->usespc                       = PETSC_FALSE;
830   ierr                               = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
831   snes->data                         = (void*)neP;
832   snes->ops->linesearchno            = SNESLineSearchNo;
833   snes->ops->linesearchnonorms       = SNESLineSearchNoNorms;
834   snes->ops->linesearchquadratic     = SNESLineSearchQuadratic_LS;
835   snes->ops->linesearchcubic         = SNESLineSearchCubic_LS;
836   snes->ops->linesearchexact         = PETSC_NULL;
837   snes->ops->linesearchtest          = PETSC_NULL;
838   snes->lsP                          = PETSC_NULL;
839   snes->ops->postcheckstep           = PETSC_NULL;
840   snes->postcheck                    = PETSC_NULL;
841   snes->ops->precheckstep            = PETSC_NULL;
842   snes->precheck                     = PETSC_NULL;
843   ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr);
844   PetscFunctionReturn(0);
845 }
846 EXTERN_C_END
847