xref: /petsc/src/snes/impls/ls/ls.c (revision b7281c8a8dad9307a4edc7c93ebc34ff865e909f)
15e76c431SBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h>
35e42470aSBarry Smith 
48a5d9ee4SBarry Smith /*
520f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
620f52e01SBarry Smith     || 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
736669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
820f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
98a5d9ee4SBarry Smith */
104a2ae208SSatish Balay #undef __FUNCT__
1104d7464bSBarry Smith #define __FUNCT__ "SNESNEWTONLSCheckLocalMin_Private"
1204d7464bSBarry Smith PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal      a1;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
16ace3abfcSBarry Smith   PetscBool      hastranspose;
178a5d9ee4SBarry Smith 
188a5d9ee4SBarry Smith   PetscFunctionBegin;
198a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2036669109SBarry Smith   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2136669109SBarry Smith   if (hastranspose) {
228a5d9ee4SBarry Smith     /* Compute || J^T F|| */
238a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
248a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
258f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
268a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2736669109SBarry Smith   } else {
2836669109SBarry Smith     Vec         work;
29ea709b57SSatish Balay     PetscScalar result;
3036669109SBarry Smith     PetscReal   wnorm;
3136669109SBarry Smith 
320298fd71SBarry Smith     ierr = VecSetRandom(W,NULL);CHKERRQ(ierr);
3336669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
376bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
3836669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
398f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
4036669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4136669109SBarry Smith   }
428a5d9ee4SBarry Smith   PetscFunctionReturn(0);
438a5d9ee4SBarry Smith }
448a5d9ee4SBarry Smith 
4574637425SBarry Smith /*
465ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4774637425SBarry Smith */
484a2ae208SSatish Balay #undef __FUNCT__
4904d7464bSBarry Smith #define __FUNCT__ "SNESNEWTONLSCheckResidual_Private"
5004d7464bSBarry Smith PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
5174637425SBarry Smith {
5274637425SBarry Smith   PetscReal      a1,a2;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
54ace3abfcSBarry Smith   PetscBool      hastranspose;
5574637425SBarry Smith 
5674637425SBarry Smith   PetscFunctionBegin;
5774637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5874637425SBarry Smith   if (hastranspose) {
5974637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6079f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6174637425SBarry Smith 
6274637425SBarry Smith     /* Compute || J^T W|| */
6374637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6474637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6574637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
6675567043SBarry Smith     if (a1 != 0.0) {
678f1a2a5eSBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
6885471664SBarry Smith     }
6974637425SBarry Smith   }
7074637425SBarry Smith   PetscFunctionReturn(0);
7174637425SBarry Smith }
7274637425SBarry Smith 
7304d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7404d965bbSLois Curfman McInnes 
7504d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7694b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7704d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
7804d965bbSLois Curfman McInnes      respectively.
7904d965bbSLois Curfman McInnes 
80fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8104d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8204d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
83fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8404d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8504d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
8604d7464bSBarry Smith      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
874b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
8804d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
8904d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9004d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9104d965bbSLois Curfman McInnes      for all nonlinear solvers.
9204d965bbSLois Curfman McInnes 
9304d965bbSLois Curfman McInnes      Another key routine is:
9404d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9504d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9604d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9704d965bbSLois Curfman McInnes      SNESSolve() if necessary.
9804d965bbSLois Curfman McInnes 
9904d965bbSLois Curfman McInnes      Additional basic routines are:
10004d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10104d965bbSLois Curfman McInnes                                       have actually been used.
10204d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
103186905e3SBarry Smith      SNESView().
10404d965bbSLois Curfman McInnes 
10504d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10604d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10704d965bbSLois Curfman McInnes      above description applies to these categories also.
10804d965bbSLois Curfman McInnes 
10904d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1105e42470aSBarry Smith /*
11104d7464bSBarry Smith    SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
11204d965bbSLois Curfman McInnes    method with a line search.
1135e76c431SBarry Smith 
11404d965bbSLois Curfman McInnes    Input Parameters:
11504d965bbSLois Curfman McInnes .  snes - the SNES context
1165e76c431SBarry Smith 
11704d965bbSLois Curfman McInnes    Output Parameter:
118c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
11904d965bbSLois Curfman McInnes 
12004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1215e76c431SBarry Smith 
1225e76c431SBarry Smith    Notes:
1235e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1245e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1255e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1265e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
127393d2d9aSLois Curfman McInnes    and Schnabel.
1285e76c431SBarry Smith */
1294a2ae208SSatish Balay #undef __FUNCT__
13004d7464bSBarry Smith #define __FUNCT__ "SNESSolve_NEWTONLS"
13104d7464bSBarry Smith PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
1325e76c431SBarry Smith {
1336849ba73SBarry Smith   PetscErrorCode      ierr;
134a7cc72afSBarry Smith   PetscInt            maxits,i,lits;
135ace3abfcSBarry Smith   PetscBool           lssucceed;
136112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
1371936c542SPeter Brune   PetscReal           fnorm,gnorm,xnorm,ynorm;
138c40d0f55SPeter Brune   Vec                 Y,X,F,G,W,FPC;
1393d4c4710SBarry Smith   KSPConvergedReason  kspreason;
1401936c542SPeter Brune   PetscBool           domainerror;
141f1c6b773SPeter Brune   SNESLineSearch      linesearch;
142c40d0f55SPeter Brune   SNESConvergedReason reason;
1435e76c431SBarry Smith 
1443a40ed3dSBarry Smith   PetscFunctionBegin;
1453d4c4710SBarry Smith   snes->numFailures            = 0;
1463d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
147184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
148184914b5SBarry Smith 
1495e42470aSBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
1505e42470aSBarry Smith   X      = snes->vec_sol;               /* solution vector */
15139e2f89bSBarry Smith   F      = snes->vec_func;              /* residual vector */
152d12b9ff3SPeter Brune   Y      = snes->vec_sol_update;        /* newton step */
153d12b9ff3SPeter Brune   G      = snes->work[0];
154d12b9ff3SPeter Brune   W      = snes->work[1];
1555e76c431SBarry Smith 
156ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1574c49b128SBarry Smith   snes->iter = 0;
15875567043SBarry Smith   snes->norm = 0.0;
159ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
1607601faf0SJed Brown   ierr       = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
161ce5a860aSPeter Brune 
162ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
163ce5a860aSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
164ce5a860aSPeter Brune     ierr = SNESApplyPC(snes,X,PETSC_NULL,PETSC_NULL,F);CHKERRQ(ierr);
165*b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
166*b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
167*b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
168*b7281c8aSPeter Brune       PetscFunctionReturn(0);
169*b7281c8aSPeter Brune     }
170*b7281c8aSPeter Brune 
171ce5a860aSPeter Brune     ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);
172ce5a860aSPeter Brune     ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
173ce5a860aSPeter Brune   } else {
174e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
17519717074SBarry Smith       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1761936c542SPeter Brune       ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1771936c542SPeter Brune       if (domainerror) {
1784936397dSBarry Smith         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1794936397dSBarry Smith         PetscFunctionReturn(0);
1804936397dSBarry Smith       }
1811aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
1821aa26658SKarl Rupp 
183e4ed7901SPeter Brune     if (!snes->norm_init_set) {
1842613ca53SBarry Smith       ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);        /* fnorm <- ||F||  */
1852613ca53SBarry Smith       ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
186189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
187189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
188189a9710SBarry Smith         PetscFunctionReturn(0);
189189a9710SBarry Smith       }
190e4ed7901SPeter Brune     } else {
191e4ed7901SPeter Brune       fnorm               = snes->norm_init;
192e4ed7901SPeter Brune       snes->norm_init_set = PETSC_FALSE;
193e4ed7901SPeter Brune     }
194ce5a860aSPeter Brune   }
195ce5a860aSPeter Brune 
196ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1975e42470aSBarry Smith   snes->norm = fnorm;
198ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
199a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
2007a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
2013f149594SLisandro Dalcin 
2023f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
2033f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
20485385478SLisandro Dalcin   /* test convergence */
20585385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
20606ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
207d034289bSBarry Smith 
2085e76c431SBarry Smith   for (i=0; i<maxits; i++) {
2095e76c431SBarry Smith 
210329e7e40SMatthew Knepley     /* Call general purpose update function */
211e7788613SBarry Smith     if (snes->ops->update) {
212e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
213de8cb200SMatthew Knepley     }
214329e7e40SMatthew Knepley 
215ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
216ce5a860aSPeter Brune     if (snes->pc) {
217ce5a860aSPeter Brune       if (snes->pcside == PC_RIGHT) {
2183024506fSPeter Brune         ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
2193024506fSPeter Brune         ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
22063e7833aSPeter Brune         ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
221c40d0f55SPeter Brune         ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr);
22263e7833aSPeter Brune         ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
223c40d0f55SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
224c40d0f55SPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
225c40d0f55SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
226c40d0f55SPeter Brune           PetscFunctionReturn(0);
227c40d0f55SPeter Brune         }
2280298fd71SBarry Smith         ierr = SNESGetFunction(snes->pc, &FPC, NULL, NULL);CHKERRQ(ierr);
229c40d0f55SPeter Brune         ierr = VecCopy(FPC, F);CHKERRQ(ierr);
230c40d0f55SPeter Brune         ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
231ce5a860aSPeter Brune       } else if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
232ce5a860aSPeter Brune         ierr = SNESApplyPC(snes,X,F,&fnorm,F);CHKERRQ(ierr);
233*b7281c8aSPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
234*b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
235*b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
236*b7281c8aSPeter Brune           PetscFunctionReturn(0);
237*b7281c8aSPeter Brune         }
238ce5a860aSPeter Brune       }
239c40d0f55SPeter Brune     }
240c40d0f55SPeter Brune 
241ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2425334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
24394b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
244d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
2453d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
2463d4c4710SBarry Smith     if (kspreason < 0) {
2473d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
248ef998cc9SBarry Smith         ierr         = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
2493d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
250ab81109fSSatish Balay         break;
2513d4c4710SBarry Smith       }
2523d4c4710SBarry Smith     }
2533d4c4710SBarry Smith     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2543f149594SLisandro Dalcin     snes->linear_its += lits;
2553f149594SLisandro Dalcin     ierr              = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
25674637425SBarry Smith 
257b0a32e0cSBarry Smith     if (PetscLogPrintInfo) {
25804d7464bSBarry Smith       ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
25985471664SBarry Smith     }
26074637425SBarry Smith 
261ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
262d12b9ff3SPeter Brune          X <- X - lambda*Y
263d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
264ea4d3ed3SLois Curfman McInnes     */
2651936c542SPeter Brune     gnorm = fnorm;
266f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
267f1c6b773SPeter Brune     ierr  = SNESLineSearchGetSuccess(linesearch, &lssucceed);CHKERRQ(ierr);
268f1c6b773SPeter Brune     ierr  = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
2691936c542SPeter Brune     ierr  = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)gnorm,(double)fnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
2704a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2711936c542SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
2721936c542SPeter Brune     if (domainerror) {
2734936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2744936397dSBarry Smith       PetscFunctionReturn(0);
2754936397dSBarry Smith     }
276a7cc72afSBarry Smith     if (!lssucceed) {
277c21ba15cSPeter Brune       if (snes->stol*xnorm > ynorm) {
278c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
279c21ba15cSPeter Brune         PetscFunctionReturn(0);
280c21ba15cSPeter Brune       }
28150ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
282ace3abfcSBarry Smith         PetscBool ismin;
283647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
28404d7464bSBarry Smith         ierr         = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2858a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2863505fcc1SBarry Smith         break;
2873505fcc1SBarry Smith       }
28850ffb88aSMatthew Knepley     }
28985385478SLisandro Dalcin     /* Monitor convergence */
290ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
29185385478SLisandro Dalcin     snes->iter = i+1;
29285385478SLisandro Dalcin     snes->norm = fnorm;
293ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
294a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
2957a03ce2fSLisandro Dalcin     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2962a997342SPeter Brune     /* Test for convergence */
297e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2983f149594SLisandro Dalcin     if (snes->reason) break;
29916c95cacSBarry Smith   }
30052392280SLois Curfman McInnes   if (i == maxits) {
301ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
30285385478SLisandro Dalcin     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
30352392280SLois Curfman McInnes   }
3043a40ed3dSBarry Smith   PetscFunctionReturn(0);
3055e76c431SBarry Smith }
30604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30704d965bbSLois Curfman McInnes /*
30804d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
30904d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
31004d965bbSLois Curfman McInnes 
31104d965bbSLois Curfman McInnes    Input Parameter:
31204d965bbSLois Curfman McInnes .  snes - the SNES context
31304d965bbSLois Curfman McInnes .  x - the solution vector
31404d965bbSLois Curfman McInnes 
31504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
31604d965bbSLois Curfman McInnes 
31704d965bbSLois Curfman McInnes    Notes:
3184b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
31904d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
32004d965bbSLois Curfman McInnes    the call to SNESSolve().
32104d965bbSLois Curfman McInnes  */
3224a2ae208SSatish Balay #undef __FUNCT__
32304d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONLS"
32404d7464bSBarry Smith PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
3255e76c431SBarry Smith {
326dfbe8321SBarry Smith   PetscErrorCode ierr;
3273a40ed3dSBarry Smith 
3283a40ed3dSBarry Smith   PetscFunctionBegin;
329fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr);
3306cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3313a40ed3dSBarry Smith   PetscFunctionReturn(0);
3325e76c431SBarry Smith }
33304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3346b8b9a38SLisandro Dalcin 
3356b8b9a38SLisandro Dalcin #undef __FUNCT__
33604d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONLS"
33704d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONLS(SNES snes)
3386b8b9a38SLisandro Dalcin {
3396b8b9a38SLisandro Dalcin   PetscFunctionBegin;
3406b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
3416b8b9a38SLisandro Dalcin }
3426b8b9a38SLisandro Dalcin 
34304d965bbSLois Curfman McInnes /*
34404d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
34504d7464bSBarry Smith    with SNESCreate_NEWTONLS().
34604d965bbSLois Curfman McInnes 
34704d965bbSLois Curfman McInnes    Input Parameter:
34804d965bbSLois Curfman McInnes .  snes - the SNES context
34904d965bbSLois Curfman McInnes 
350de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
35104d965bbSLois Curfman McInnes  */
3524a2ae208SSatish Balay #undef __FUNCT__
35304d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONLS"
35404d7464bSBarry Smith PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
3555e76c431SBarry Smith {
356dfbe8321SBarry Smith   PetscErrorCode ierr;
3573a40ed3dSBarry Smith 
3583a40ed3dSBarry Smith   PetscFunctionBegin;
35904d7464bSBarry Smith   ierr = SNESReset_NEWTONLS(snes);CHKERRQ(ierr);
3605c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3613a40ed3dSBarry Smith   PetscFunctionReturn(0);
3625e76c431SBarry Smith }
36304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
36494298653SBarry Smith 
365329e7e40SMatthew Knepley /*
36604d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
36704d965bbSLois Curfman McInnes 
36804d965bbSLois Curfman McInnes    Input Parameters:
36904d965bbSLois Curfman McInnes .  SNES - the SNES context
37004d965bbSLois Curfman McInnes .  viewer - visualization context
37104d965bbSLois Curfman McInnes 
37204d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
37304d965bbSLois Curfman McInnes */
3744a2ae208SSatish Balay #undef __FUNCT__
37504d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONLS"
37604d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer)
377a935fc98SLois Curfman McInnes {
378dfbe8321SBarry Smith   PetscErrorCode ierr;
379ace3abfcSBarry Smith   PetscBool      iascii;
380a935fc98SLois Curfman McInnes 
3813a40ed3dSBarry Smith   PetscFunctionBegin;
382251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
38332077d6dSBarry Smith   if (iascii) {
38450f6de3fSJed Brown   }
38550f6de3fSJed Brown   PetscFunctionReturn(0);
38650f6de3fSJed Brown }
38750f6de3fSJed Brown 
38804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
38904d965bbSLois Curfman McInnes /*
39004d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
39104d965bbSLois Curfman McInnes 
39204d965bbSLois Curfman McInnes    Input Parameter:
39304d965bbSLois Curfman McInnes .  snes - the SNES context
39404d965bbSLois Curfman McInnes 
39504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
39604d965bbSLois Curfman McInnes */
3974a2ae208SSatish Balay #undef __FUNCT__
39804d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONLS"
39904d7464bSBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes)
4005e42470aSBarry Smith {
401dfbe8321SBarry Smith   PetscErrorCode ierr;
402f1c6b773SPeter Brune   SNESLineSearch linesearch;
4035e42470aSBarry Smith 
4043a40ed3dSBarry Smith   PetscFunctionBegin;
40504d7464bSBarry Smith   ierr = PetscOptionsHead("SNESNEWTONLS options");CHKERRQ(ierr);
406b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4079e764e56SPeter Brune   /* set the default line search type */
4089e764e56SPeter Brune   if (!snes->linesearch) {
4097601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
4101a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
4119e764e56SPeter Brune   }
4123a40ed3dSBarry Smith   PetscFunctionReturn(0);
4135e42470aSBarry Smith }
4144827ddcaSBarry Smith 
41504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
416ebe3b25bSBarry Smith /*MC
41704d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
41804d965bbSLois Curfman McInnes 
419ebe3b25bSBarry Smith    Options Database:
4205a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
4215a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
4225a9b6599SPeter Brune .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
42382d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
42482d26c24SPeter Brune .   -snes_linesearch_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)
4255a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
42682d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
42782d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
428acbee50cSBarry Smith 
429ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
43004d965bbSLois Curfman McInnes 
431ee3001cbSBarry Smith    Level: beginner
432ee3001cbSBarry Smith 
43304d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder()
4345a9b6599SPeter Brune            SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms()
435ebe3b25bSBarry Smith 
436ebe3b25bSBarry Smith M*/
4374a2ae208SSatish Balay #undef __FUNCT__
43804d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONLS"
4398cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
4405e42470aSBarry Smith {
441dfbe8321SBarry Smith   PetscErrorCode ierr;
44204d7464bSBarry Smith   SNES_NEWTONLS  *neP;
4435e42470aSBarry Smith 
4443a40ed3dSBarry Smith   PetscFunctionBegin;
44504d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
44604d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
44704d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
44804d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
44904d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
45004d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
4515e42470aSBarry Smith 
45242f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
4536b6c8463SPeter Brune   snes->usespc  = PETSC_TRUE;
45404d7464bSBarry Smith   ierr          = PetscNewLog(snes,SNES_NEWTONLS,&neP);CHKERRQ(ierr);
4555e42470aSBarry Smith   snes->data    = (void*)neP;
4563a40ed3dSBarry Smith   PetscFunctionReturn(0);
4575e42470aSBarry Smith }
458