xref: /petsc/src/snes/impls/ls/ls.c (revision f6c5ca7819faa18d8a3a15fb5c97bf2c315aa240)
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"
12*f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,PetscReal fnorm,PetscBool  *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal      a1;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
16ace3abfcSBarry Smith   PetscBool      hastranspose;
17*f6c5ca78SBarry Smith   Vec            W;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2136669109SBarry Smith   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
22*f6c5ca78SBarry Smith   ierr   = VecDuplicate(F,&W);CHKERRQ(ierr);
2336669109SBarry Smith   if (hastranspose) {
248a5d9ee4SBarry Smith     /* Compute || J^T F|| */
258a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
268a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
278f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
288a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2936669109SBarry Smith   } else {
3036669109SBarry Smith     Vec         work;
31ea709b57SSatish Balay     PetscScalar result;
3236669109SBarry Smith     PetscReal   wnorm;
3336669109SBarry Smith 
340298fd71SBarry Smith     ierr = VecSetRandom(W,NULL);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3736669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3836669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
396bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
4036669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
418f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
4236669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4336669109SBarry Smith   }
44*f6c5ca78SBarry Smith   ierr = VecDestroy(&W);CHKERRQ(ierr);
458a5d9ee4SBarry Smith   PetscFunctionReturn(0);
468a5d9ee4SBarry Smith }
478a5d9ee4SBarry Smith 
4874637425SBarry Smith /*
495ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
5074637425SBarry Smith */
514a2ae208SSatish Balay #undef __FUNCT__
5204d7464bSBarry Smith #define __FUNCT__ "SNESNEWTONLSCheckResidual_Private"
53*f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X)
5474637425SBarry Smith {
5574637425SBarry Smith   PetscReal      a1,a2;
56dfbe8321SBarry Smith   PetscErrorCode ierr;
57ace3abfcSBarry Smith   PetscBool      hastranspose;
5874637425SBarry Smith 
5974637425SBarry Smith   PetscFunctionBegin;
6074637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
6174637425SBarry Smith   if (hastranspose) {
62*f6c5ca78SBarry Smith     Vec   W1,W2;
63*f6c5ca78SBarry Smith 
64*f6c5ca78SBarry Smith     ierr = VecDuplicate(F,&W1);CHKERRQ(ierr);
65*f6c5ca78SBarry Smith     ierr = VecDuplicate(F,&W2);CHKERRQ(ierr);
6674637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6779f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6874637425SBarry Smith 
6974637425SBarry Smith     /* Compute || J^T W|| */
7074637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
7174637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
7274637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
7375567043SBarry Smith     if (a1 != 0.0) {
748f1a2a5eSBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
7585471664SBarry Smith     }
76*f6c5ca78SBarry Smith     ierr = VecDestroy(&W1);CHKERRQ(ierr);
77*f6c5ca78SBarry Smith     ierr = VecDestroy(&W2);CHKERRQ(ierr);
7874637425SBarry Smith   }
7974637425SBarry Smith   PetscFunctionReturn(0);
8074637425SBarry Smith }
8174637425SBarry Smith 
8204d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
8304d965bbSLois Curfman McInnes 
8404d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
8594b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
8604d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8704d965bbSLois Curfman McInnes      respectively.
8804d965bbSLois Curfman McInnes 
89fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
9004d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
9104d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
92fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
9304d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
9404d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
9504d7464bSBarry Smith      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
964b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
9704d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9804d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9904d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
10004d965bbSLois Curfman McInnes      for all nonlinear solvers.
10104d965bbSLois Curfman McInnes 
10204d965bbSLois Curfman McInnes      Another key routine is:
10304d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
10404d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
10504d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
10604d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10704d965bbSLois Curfman McInnes 
10804d965bbSLois Curfman McInnes      Additional basic routines are:
10904d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
11004d965bbSLois Curfman McInnes                                       have actually been used.
11104d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
112186905e3SBarry Smith      SNESView().
11304d965bbSLois Curfman McInnes 
11404d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
11504d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
11604d965bbSLois Curfman McInnes      above description applies to these categories also.
11704d965bbSLois Curfman McInnes 
11804d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1195e42470aSBarry Smith /*
12004d7464bSBarry Smith    SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
12104d965bbSLois Curfman McInnes    method with a line search.
1225e76c431SBarry Smith 
12304d965bbSLois Curfman McInnes    Input Parameters:
12404d965bbSLois Curfman McInnes .  snes - the SNES context
1255e76c431SBarry Smith 
12604d965bbSLois Curfman McInnes    Output Parameter:
127c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12804d965bbSLois Curfman McInnes 
12904d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1305e76c431SBarry Smith 
1315e76c431SBarry Smith    Notes:
1325e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1335e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1345e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1355e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
136393d2d9aSLois Curfman McInnes    and Schnabel.
1375e76c431SBarry Smith */
1384a2ae208SSatish Balay #undef __FUNCT__
13904d7464bSBarry Smith #define __FUNCT__ "SNESSolve_NEWTONLS"
14004d7464bSBarry Smith PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
1415e76c431SBarry Smith {
1426849ba73SBarry Smith   PetscErrorCode      ierr;
143a7cc72afSBarry Smith   PetscInt            maxits,i,lits;
144ace3abfcSBarry Smith   PetscBool           lssucceed;
1451936c542SPeter Brune   PetscReal           fnorm,gnorm,xnorm,ynorm;
146*f6c5ca78SBarry Smith   Vec                 Y,X,F;
1473d4c4710SBarry Smith   KSPConvergedReason  kspreason;
1481936c542SPeter Brune   PetscBool           domainerror;
149f1c6b773SPeter Brune   SNESLineSearch      linesearch;
150c40d0f55SPeter Brune   SNESConvergedReason reason;
1515e76c431SBarry Smith 
1523a40ed3dSBarry Smith   PetscFunctionBegin;
1533d4c4710SBarry Smith   snes->numFailures            = 0;
1543d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
155184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
156184914b5SBarry Smith 
1575e42470aSBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
1585e42470aSBarry Smith   X      = snes->vec_sol;               /* solution vector */
15939e2f89bSBarry Smith   F      = snes->vec_func;              /* residual vector */
160d12b9ff3SPeter Brune   Y      = snes->vec_sol_update;        /* newton step */
1615e76c431SBarry Smith 
162e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1634c49b128SBarry Smith   snes->iter = 0;
16475567043SBarry Smith   snes->norm = 0.0;
165e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
1667601faf0SJed Brown   ierr       = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
167ce5a860aSPeter Brune 
168ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
169ce5a860aSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
170be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
171b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
172b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
173b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
174b7281c8aSPeter Brune       PetscFunctionReturn(0);
175b7281c8aSPeter Brune     }
176b7281c8aSPeter Brune 
177ce5a860aSPeter Brune     ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);
178ce5a860aSPeter Brune     ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
179ce5a860aSPeter Brune   } else {
180e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
18119717074SBarry Smith       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1821936c542SPeter Brune       ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1831936c542SPeter Brune       if (domainerror) {
1844936397dSBarry Smith         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1854936397dSBarry Smith         PetscFunctionReturn(0);
1864936397dSBarry Smith       }
1871aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
188c1c75074SPeter Brune   }
1891aa26658SKarl Rupp 
190c1c75074SPeter Brune   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);        /* fnorm <- ||F||  */
191189a9710SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) {
192189a9710SBarry Smith     snes->reason = SNES_DIVERGED_FNORM_NAN;
193189a9710SBarry Smith     PetscFunctionReturn(0);
194189a9710SBarry Smith   }
195ce5a860aSPeter Brune 
196e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1975e42470aSBarry Smith   snes->norm = fnorm;
198e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
199a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
2007a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
2013f149594SLisandro Dalcin 
20285385478SLisandro Dalcin   /* test convergence */
20385385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
20406ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
205d034289bSBarry Smith 
2065e76c431SBarry Smith   for (i=0; i<maxits; i++) {
2075e76c431SBarry Smith 
208329e7e40SMatthew Knepley     /* Call general purpose update function */
209e7788613SBarry Smith     if (snes->ops->update) {
210e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
211de8cb200SMatthew Knepley     }
212329e7e40SMatthew Knepley 
213ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
214ce5a860aSPeter Brune     if (snes->pc) {
215ce5a860aSPeter Brune       if (snes->pcside == PC_RIGHT) {
2163024506fSPeter Brune         ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
21763e7833aSPeter Brune         ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
218c40d0f55SPeter Brune         ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr);
21963e7833aSPeter Brune         ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
220c40d0f55SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
221c40d0f55SPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
222c40d0f55SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
223c40d0f55SPeter Brune           PetscFunctionReturn(0);
224c40d0f55SPeter Brune         }
225be95d8f1SBarry Smith         ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
226ce5a860aSPeter Brune       } else if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
227be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,F,F);CHKERRQ(ierr);
228b7281c8aSPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
229b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
230b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
231b7281c8aSPeter Brune           PetscFunctionReturn(0);
232b7281c8aSPeter Brune         }
233ce5a860aSPeter Brune       }
234c40d0f55SPeter Brune     }
235c40d0f55SPeter Brune 
236ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
237d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
23823ee1639SBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
239d4211eb9SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
2403d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
2413d4c4710SBarry Smith     if (kspreason < 0) {
2423d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
243ef998cc9SBarry 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);
2443d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
245ab81109fSSatish Balay         break;
2463d4c4710SBarry Smith       }
2473d4c4710SBarry Smith     }
2483d4c4710SBarry Smith     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2493f149594SLisandro Dalcin     snes->linear_its += lits;
2503f149594SLisandro Dalcin     ierr              = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
25174637425SBarry Smith 
252b0a32e0cSBarry Smith     if (PetscLogPrintInfo) {
253*f6c5ca78SBarry Smith       ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y);CHKERRQ(ierr);
25485471664SBarry Smith     }
25574637425SBarry Smith 
256ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
257d12b9ff3SPeter Brune          X <- X - lambda*Y
258d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
259ea4d3ed3SLois Curfman McInnes     */
2601936c542SPeter Brune     gnorm = fnorm;
261f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
262f1c6b773SPeter Brune     ierr  = SNESLineSearchGetSuccess(linesearch, &lssucceed);CHKERRQ(ierr);
263f1c6b773SPeter Brune     ierr  = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
2641936c542SPeter 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);
2654a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2661936c542SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
2671936c542SPeter Brune     if (domainerror) {
2684936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2694936397dSBarry Smith       PetscFunctionReturn(0);
2704936397dSBarry Smith     }
271a7cc72afSBarry Smith     if (!lssucceed) {
272c21ba15cSPeter Brune       if (snes->stol*xnorm > ynorm) {
273c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
274c21ba15cSPeter Brune         PetscFunctionReturn(0);
275c21ba15cSPeter Brune       }
27650ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
277ace3abfcSBarry Smith         PetscBool ismin;
278647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
279*f6c5ca78SBarry Smith         ierr         = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,fnorm,&ismin);CHKERRQ(ierr);
2808a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2813505fcc1SBarry Smith         break;
2823505fcc1SBarry Smith       }
28350ffb88aSMatthew Knepley     }
28485385478SLisandro Dalcin     /* Monitor convergence */
285e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
28685385478SLisandro Dalcin     snes->iter = i+1;
28785385478SLisandro Dalcin     snes->norm = fnorm;
288e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
289a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
2907a03ce2fSLisandro Dalcin     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2912a997342SPeter Brune     /* Test for convergence */
292e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2933f149594SLisandro Dalcin     if (snes->reason) break;
29416c95cacSBarry Smith   }
29552392280SLois Curfman McInnes   if (i == maxits) {
296ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
29785385478SLisandro Dalcin     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
29852392280SLois Curfman McInnes   }
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
3005e76c431SBarry Smith }
30104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30204d965bbSLois Curfman McInnes /*
30304d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
30404d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
30504d965bbSLois Curfman McInnes 
30604d965bbSLois Curfman McInnes    Input Parameter:
30704d965bbSLois Curfman McInnes .  snes - the SNES context
30804d965bbSLois Curfman McInnes .  x - the solution vector
30904d965bbSLois Curfman McInnes 
31004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
31104d965bbSLois Curfman McInnes 
31204d965bbSLois Curfman McInnes    Notes:
3134b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
31404d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
31504d965bbSLois Curfman McInnes    the call to SNESSolve().
31604d965bbSLois Curfman McInnes  */
3174a2ae208SSatish Balay #undef __FUNCT__
31804d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONLS"
31904d7464bSBarry Smith PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
3205e76c431SBarry Smith {
321dfbe8321SBarry Smith   PetscErrorCode ierr;
3223a40ed3dSBarry Smith 
3233a40ed3dSBarry Smith   PetscFunctionBegin;
3246cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3256c67d002SPeter Brune   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
3275e76c431SBarry Smith }
32804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3296b8b9a38SLisandro Dalcin 
3306b8b9a38SLisandro Dalcin #undef __FUNCT__
33104d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONLS"
33204d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONLS(SNES snes)
3336b8b9a38SLisandro Dalcin {
3346b8b9a38SLisandro Dalcin   PetscFunctionBegin;
3356b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
3366b8b9a38SLisandro Dalcin }
3376b8b9a38SLisandro Dalcin 
33804d965bbSLois Curfman McInnes /*
33904d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
34004d7464bSBarry Smith    with SNESCreate_NEWTONLS().
34104d965bbSLois Curfman McInnes 
34204d965bbSLois Curfman McInnes    Input Parameter:
34304d965bbSLois Curfman McInnes .  snes - the SNES context
34404d965bbSLois Curfman McInnes 
345de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
34604d965bbSLois Curfman McInnes  */
3474a2ae208SSatish Balay #undef __FUNCT__
34804d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONLS"
34904d7464bSBarry Smith PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
3505e76c431SBarry Smith {
351dfbe8321SBarry Smith   PetscErrorCode ierr;
3523a40ed3dSBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
35404d7464bSBarry Smith   ierr = SNESReset_NEWTONLS(snes);CHKERRQ(ierr);
3555c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3563a40ed3dSBarry Smith   PetscFunctionReturn(0);
3575e76c431SBarry Smith }
35804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
35994298653SBarry Smith 
360329e7e40SMatthew Knepley /*
36104d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
36204d965bbSLois Curfman McInnes 
36304d965bbSLois Curfman McInnes    Input Parameters:
36404d965bbSLois Curfman McInnes .  SNES - the SNES context
36504d965bbSLois Curfman McInnes .  viewer - visualization context
36604d965bbSLois Curfman McInnes 
36704d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
36804d965bbSLois Curfman McInnes */
3694a2ae208SSatish Balay #undef __FUNCT__
37004d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONLS"
37104d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer)
372a935fc98SLois Curfman McInnes {
373dfbe8321SBarry Smith   PetscErrorCode ierr;
374ace3abfcSBarry Smith   PetscBool      iascii;
375a935fc98SLois Curfman McInnes 
3763a40ed3dSBarry Smith   PetscFunctionBegin;
377251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
37832077d6dSBarry Smith   if (iascii) {
37950f6de3fSJed Brown   }
38050f6de3fSJed Brown   PetscFunctionReturn(0);
38150f6de3fSJed Brown }
38250f6de3fSJed Brown 
38304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
38404d965bbSLois Curfman McInnes /*
38504d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
38604d965bbSLois Curfman McInnes 
38704d965bbSLois Curfman McInnes    Input Parameter:
38804d965bbSLois Curfman McInnes .  snes - the SNES context
38904d965bbSLois Curfman McInnes 
39004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
39104d965bbSLois Curfman McInnes */
3924a2ae208SSatish Balay #undef __FUNCT__
39304d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONLS"
3948c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(PetscOptions *PetscOptionsObject,SNES snes)
3955e42470aSBarry Smith {
396dfbe8321SBarry Smith   PetscErrorCode ierr;
397f1c6b773SPeter Brune   SNESLineSearch linesearch;
3985e42470aSBarry Smith 
3993a40ed3dSBarry Smith   PetscFunctionBegin;
4009e764e56SPeter Brune   if (!snes->linesearch) {
4017601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
4021a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
4039e764e56SPeter Brune   }
4043a40ed3dSBarry Smith   PetscFunctionReturn(0);
4055e42470aSBarry Smith }
4064827ddcaSBarry Smith 
40704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
408ebe3b25bSBarry Smith /*MC
40904d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
41004d965bbSLois Curfman McInnes 
411ebe3b25bSBarry Smith    Options Database:
4125a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
4135a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
4145a9b6599SPeter Brune .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
41582d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
41682d26c24SPeter 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)
4175a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
41882d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
41982d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
420acbee50cSBarry Smith 
421ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
42204d965bbSLois Curfman McInnes 
423ee3001cbSBarry Smith    Level: beginner
424ee3001cbSBarry Smith 
42504d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder()
4265a9b6599SPeter Brune            SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms()
427ebe3b25bSBarry Smith 
428ebe3b25bSBarry Smith M*/
4294a2ae208SSatish Balay #undef __FUNCT__
43004d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONLS"
4318cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
4325e42470aSBarry Smith {
433dfbe8321SBarry Smith   PetscErrorCode ierr;
43404d7464bSBarry Smith   SNES_NEWTONLS  *neP;
4355e42470aSBarry Smith 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
43704d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
43804d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
43904d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
44004d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
44104d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
44204d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
4435e42470aSBarry Smith 
4446c67d002SPeter Brune   snes->pcside  = PC_RIGHT;
44542f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
4466b6c8463SPeter Brune   snes->usespc  = PETSC_TRUE;
447b00a9115SJed Brown   ierr          = PetscNewLog(snes,&neP);CHKERRQ(ierr);
4485e42470aSBarry Smith   snes->data    = (void*)neP;
4493a40ed3dSBarry Smith   PetscFunctionReturn(0);
4505e42470aSBarry Smith }
451