xref: /petsc/src/snes/impls/ls/ls.c (revision c579b30025f7481884435f5d6baf5e954044f830)
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"
12f6c5ca78SBarry 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;
17f6c5ca78SBarry Smith   Vec            W;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2136669109SBarry Smith   ierr   = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
22f6c5ca78SBarry 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   }
44f6c5ca78SBarry 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"
53f6c5ca78SBarry 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) {
62f6c5ca78SBarry Smith     Vec   W1,W2;
63f6c5ca78SBarry Smith 
64f6c5ca78SBarry Smith     ierr = VecDuplicate(F,&W1);CHKERRQ(ierr);
65f6c5ca78SBarry 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     }
76f6c5ca78SBarry Smith     ierr = VecDestroy(&W1);CHKERRQ(ierr);
77f6c5ca78SBarry 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;
146f6c5ca78SBarry 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;
153*c579b300SPatrick Farrell 
154*c579b300SPatrick Farrell   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
155*c579b300SPatrick Farrell     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
156*c579b300SPatrick Farrell   }
157*c579b300SPatrick Farrell 
1583d4c4710SBarry Smith   snes->numFailures            = 0;
1593d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
160184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
161184914b5SBarry Smith 
1625e42470aSBarry Smith   maxits = snes->max_its;               /* maximum number of iterations */
1635e42470aSBarry Smith   X      = snes->vec_sol;               /* solution vector */
16439e2f89bSBarry Smith   F      = snes->vec_func;              /* residual vector */
165d12b9ff3SPeter Brune   Y      = snes->vec_sol_update;        /* newton step */
1665e76c431SBarry Smith 
167e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
1684c49b128SBarry Smith   snes->iter = 0;
16975567043SBarry Smith   snes->norm = 0.0;
170e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
1717601faf0SJed Brown   ierr       = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
172ce5a860aSPeter Brune 
173ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
174ce5a860aSPeter Brune   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
175be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
176b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
177b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
178b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
179b7281c8aSPeter Brune       PetscFunctionReturn(0);
180b7281c8aSPeter Brune     }
181b7281c8aSPeter Brune 
182ce5a860aSPeter Brune     ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);
183ce5a860aSPeter Brune     ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
184ce5a860aSPeter Brune   } else {
185e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
18619717074SBarry Smith       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1871936c542SPeter Brune       ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1881936c542SPeter Brune       if (domainerror) {
1894936397dSBarry Smith         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1904936397dSBarry Smith         PetscFunctionReturn(0);
1914936397dSBarry Smith       }
1921aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
193c1c75074SPeter Brune   }
1941aa26658SKarl Rupp 
195c1c75074SPeter Brune   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);        /* fnorm <- ||F||  */
196189a9710SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) {
197189a9710SBarry Smith     snes->reason = SNES_DIVERGED_FNORM_NAN;
198189a9710SBarry Smith     PetscFunctionReturn(0);
199189a9710SBarry Smith   }
200ce5a860aSPeter Brune 
201e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
2025e42470aSBarry Smith   snes->norm = fnorm;
203e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
204a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
2057a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
2063f149594SLisandro Dalcin 
20785385478SLisandro Dalcin   /* test convergence */
20885385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
20906ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
210d034289bSBarry Smith 
2115e76c431SBarry Smith   for (i=0; i<maxits; i++) {
2125e76c431SBarry Smith 
213329e7e40SMatthew Knepley     /* Call general purpose update function */
214e7788613SBarry Smith     if (snes->ops->update) {
215e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
216de8cb200SMatthew Knepley     }
217329e7e40SMatthew Knepley 
218ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
219ce5a860aSPeter Brune     if (snes->pc) {
220ce5a860aSPeter Brune       if (snes->pcside == PC_RIGHT) {
2213024506fSPeter Brune         ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
22263e7833aSPeter Brune         ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
223c40d0f55SPeter Brune         ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr);
22463e7833aSPeter Brune         ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
225c40d0f55SPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
226c40d0f55SPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
227c40d0f55SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
228c40d0f55SPeter Brune           PetscFunctionReturn(0);
229c40d0f55SPeter Brune         }
230be95d8f1SBarry Smith         ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr);
231ce5a860aSPeter Brune       } else if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
232be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,F,F);CHKERRQ(ierr);
233b7281c8aSPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
234b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
235b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
236b7281c8aSPeter Brune           PetscFunctionReturn(0);
237b7281c8aSPeter Brune         }
238ce5a860aSPeter Brune       }
239c40d0f55SPeter Brune     }
240c40d0f55SPeter Brune 
241ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
242d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
24323ee1639SBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);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) {
258f6c5ca78SBarry Smith       ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y);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;
284f6c5ca78SBarry Smith         ierr         = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,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 */
290e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
29185385478SLisandro Dalcin     snes->iter = i+1;
29285385478SLisandro Dalcin     snes->norm = fnorm;
293e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((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;
3296cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3306c67d002SPeter Brune   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
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"
3998c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(PetscOptions *PetscOptionsObject,SNES snes)
4005e42470aSBarry Smith {
401dfbe8321SBarry Smith   PetscErrorCode ierr;
402f1c6b773SPeter Brune   SNESLineSearch linesearch;
4035e42470aSBarry Smith 
4043a40ed3dSBarry Smith   PetscFunctionBegin;
4059e764e56SPeter Brune   if (!snes->linesearch) {
4067601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
4071a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
4089e764e56SPeter Brune   }
4093a40ed3dSBarry Smith   PetscFunctionReturn(0);
4105e42470aSBarry Smith }
4114827ddcaSBarry Smith 
41204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
413ebe3b25bSBarry Smith /*MC
41404d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
41504d965bbSLois Curfman McInnes 
416ebe3b25bSBarry Smith    Options Database:
4175a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
4185a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
4195a9b6599SPeter Brune .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
42082d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
42182d26c24SPeter 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)
4225a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
42382d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
42482d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
425acbee50cSBarry Smith 
426ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
42704d965bbSLois Curfman McInnes 
428ee3001cbSBarry Smith    Level: beginner
429ee3001cbSBarry Smith 
43004d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder()
4315a9b6599SPeter Brune            SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms()
432ebe3b25bSBarry Smith 
433ebe3b25bSBarry Smith M*/
4344a2ae208SSatish Balay #undef __FUNCT__
43504d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONLS"
4368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
4375e42470aSBarry Smith {
438dfbe8321SBarry Smith   PetscErrorCode ierr;
43904d7464bSBarry Smith   SNES_NEWTONLS  *neP;
4405e42470aSBarry Smith 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
44204d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
44304d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
44404d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
44504d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
44604d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
44704d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
4485e42470aSBarry Smith 
4496c67d002SPeter Brune   snes->pcside  = PC_RIGHT;
45042f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
4516b6c8463SPeter Brune   snes->usespc  = PETSC_TRUE;
452b00a9115SJed Brown   ierr          = PetscNewLog(snes,&neP);CHKERRQ(ierr);
4535e42470aSBarry Smith   snes->data    = (void*)neP;
4543a40ed3dSBarry Smith   PetscFunctionReturn(0);
4555e42470aSBarry Smith }
456