xref: /petsc/src/snes/impls/ls/ls.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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 
1564c49b128SBarry Smith   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1574c49b128SBarry Smith   snes->iter = 0;
15875567043SBarry Smith   snes->norm = 0.0;
1594c49b128SBarry Smith   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
160e4ed7901SPeter Brune   ierr       = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
161e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
16219717074SBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1631936c542SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1641936c542SPeter Brune     if (domainerror) {
1654936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1664936397dSBarry Smith       PetscFunctionReturn(0);
1674936397dSBarry Smith     }
1681aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
1691aa26658SKarl Rupp 
170e4ed7901SPeter Brune   if (!snes->norm_init_set) {
1712613ca53SBarry Smith     ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);        /* fnorm <- ||F||  */
1722613ca53SBarry Smith     ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
173189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
174189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
175189a9710SBarry Smith       PetscFunctionReturn(0);
176189a9710SBarry Smith     }
177e4ed7901SPeter Brune   } else {
178e4ed7901SPeter Brune     fnorm               = snes->norm_init;
179e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
180e4ed7901SPeter Brune   }
1810f5bd95cSBarry Smith   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1825e42470aSBarry Smith   snes->norm = fnorm;
1830f5bd95cSBarry Smith   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
184a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
1857a03ce2fSLisandro Dalcin   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1863f149594SLisandro Dalcin 
1873f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
1883f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
18985385478SLisandro Dalcin   /* test convergence */
19085385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
19106ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
192d034289bSBarry Smith 
1935e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1945e76c431SBarry Smith 
195329e7e40SMatthew Knepley     /* Call general purpose update function */
196e7788613SBarry Smith     if (snes->ops->update) {
197e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
198de8cb200SMatthew Knepley     }
199329e7e40SMatthew Knepley 
200c40d0f55SPeter Brune     /* apply the nonlinear preconditioner if it's right preconditioned */
201c40d0f55SPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
2023024506fSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
2033024506fSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
20463e7833aSPeter Brune       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
205c40d0f55SPeter Brune       ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr);
20663e7833aSPeter Brune       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr);
207c40d0f55SPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
208c40d0f55SPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
209c40d0f55SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
210c40d0f55SPeter Brune         PetscFunctionReturn(0);
211c40d0f55SPeter Brune       }
2120298fd71SBarry Smith       ierr = SNESGetFunction(snes->pc, &FPC, NULL, NULL);CHKERRQ(ierr);
213c40d0f55SPeter Brune       ierr = VecCopy(FPC, F);CHKERRQ(ierr);
214c40d0f55SPeter Brune       ierr = SNESGetFunctionNorm(snes->pc, &fnorm);CHKERRQ(ierr);
215c40d0f55SPeter Brune     }
216c40d0f55SPeter Brune 
217ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2185334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
21994b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
22071f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
2213d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
2223d4c4710SBarry Smith     if (kspreason < 0) {
2233d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
224ef998cc9SBarry 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);
2253d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
226ab81109fSSatish Balay         break;
2273d4c4710SBarry Smith       }
2283d4c4710SBarry Smith     }
2293d4c4710SBarry Smith     ierr              = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2303f149594SLisandro Dalcin     snes->linear_its += lits;
2313f149594SLisandro Dalcin     ierr              = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
23274637425SBarry Smith 
233b0a32e0cSBarry Smith     if (PetscLogPrintInfo) {
23404d7464bSBarry Smith       ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
23585471664SBarry Smith     }
23674637425SBarry Smith 
237ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
238d12b9ff3SPeter Brune          X <- X - lambda*Y
239d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
240ea4d3ed3SLois Curfman McInnes     */
2411936c542SPeter Brune     gnorm = fnorm;
242f1c6b773SPeter Brune     ierr  = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
243f1c6b773SPeter Brune     ierr  = SNESLineSearchGetSuccess(linesearch, &lssucceed);CHKERRQ(ierr);
244f1c6b773SPeter Brune     ierr  = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
2451936c542SPeter 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);
2464a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2471936c542SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
2481936c542SPeter Brune     if (domainerror) {
2494936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2504936397dSBarry Smith       PetscFunctionReturn(0);
2514936397dSBarry Smith     }
252a7cc72afSBarry Smith     if (!lssucceed) {
253c21ba15cSPeter Brune       if (snes->stol*xnorm > ynorm) {
254c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
255c21ba15cSPeter Brune         PetscFunctionReturn(0);
256c21ba15cSPeter Brune       }
25750ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
258ace3abfcSBarry Smith         PetscBool ismin;
259647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
26004d7464bSBarry Smith         ierr         = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2618a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2623505fcc1SBarry Smith         break;
2633505fcc1SBarry Smith       }
26450ffb88aSMatthew Knepley     }
26585385478SLisandro Dalcin     /* Monitor convergence */
26685385478SLisandro Dalcin     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
26785385478SLisandro Dalcin     snes->iter = i+1;
26885385478SLisandro Dalcin     snes->norm = fnorm;
26985385478SLisandro Dalcin     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
270a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
2717a03ce2fSLisandro Dalcin     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
2722a997342SPeter Brune     /* Test for convergence */
273e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2743f149594SLisandro Dalcin     if (snes->reason) break;
27516c95cacSBarry Smith   }
27652392280SLois Curfman McInnes   if (i == maxits) {
277ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
27885385478SLisandro Dalcin     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
27952392280SLois Curfman McInnes   }
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
2815e76c431SBarry Smith }
28204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28304d965bbSLois Curfman McInnes /*
28404d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
28504d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
28604d965bbSLois Curfman McInnes 
28704d965bbSLois Curfman McInnes    Input Parameter:
28804d965bbSLois Curfman McInnes .  snes - the SNES context
28904d965bbSLois Curfman McInnes .  x - the solution vector
29004d965bbSLois Curfman McInnes 
29104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
29204d965bbSLois Curfman McInnes 
29304d965bbSLois Curfman McInnes    Notes:
2944b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
29504d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
29604d965bbSLois Curfman McInnes    the call to SNESSolve().
29704d965bbSLois Curfman McInnes  */
2984a2ae208SSatish Balay #undef __FUNCT__
29904d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONLS"
30004d7464bSBarry Smith PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
3015e76c431SBarry Smith {
302dfbe8321SBarry Smith   PetscErrorCode ierr;
3033a40ed3dSBarry Smith 
3043a40ed3dSBarry Smith   PetscFunctionBegin;
305d12b9ff3SPeter Brune   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
3066cab3a1bSJed Brown   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3073a40ed3dSBarry Smith   PetscFunctionReturn(0);
3085e76c431SBarry Smith }
30904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3106b8b9a38SLisandro Dalcin 
3116b8b9a38SLisandro Dalcin #undef __FUNCT__
31204d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONLS"
31304d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONLS(SNES snes)
3146b8b9a38SLisandro Dalcin {
3156b8b9a38SLisandro Dalcin   PetscFunctionBegin;
3166b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
3176b8b9a38SLisandro Dalcin }
3186b8b9a38SLisandro Dalcin 
31904d965bbSLois Curfman McInnes /*
32004d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
32104d7464bSBarry Smith    with SNESCreate_NEWTONLS().
32204d965bbSLois Curfman McInnes 
32304d965bbSLois Curfman McInnes    Input Parameter:
32404d965bbSLois Curfman McInnes .  snes - the SNES context
32504d965bbSLois Curfman McInnes 
326de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
32704d965bbSLois Curfman McInnes  */
3284a2ae208SSatish Balay #undef __FUNCT__
32904d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONLS"
33004d7464bSBarry Smith PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
3315e76c431SBarry Smith {
332dfbe8321SBarry Smith   PetscErrorCode ierr;
3333a40ed3dSBarry Smith 
3343a40ed3dSBarry Smith   PetscFunctionBegin;
33504d7464bSBarry Smith   ierr = SNESReset_NEWTONLS(snes);CHKERRQ(ierr);
3365c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3373a40ed3dSBarry Smith   PetscFunctionReturn(0);
3385e76c431SBarry Smith }
33904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
34094298653SBarry Smith 
341329e7e40SMatthew Knepley /*
34204d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
34304d965bbSLois Curfman McInnes 
34404d965bbSLois Curfman McInnes    Input Parameters:
34504d965bbSLois Curfman McInnes .  SNES - the SNES context
34604d965bbSLois Curfman McInnes .  viewer - visualization context
34704d965bbSLois Curfman McInnes 
34804d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
34904d965bbSLois Curfman McInnes */
3504a2ae208SSatish Balay #undef __FUNCT__
35104d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONLS"
35204d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer)
353a935fc98SLois Curfman McInnes {
354dfbe8321SBarry Smith   PetscErrorCode ierr;
355ace3abfcSBarry Smith   PetscBool      iascii;
356a935fc98SLois Curfman McInnes 
3573a40ed3dSBarry Smith   PetscFunctionBegin;
358251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
35932077d6dSBarry Smith   if (iascii) {
36050f6de3fSJed Brown   }
36150f6de3fSJed Brown   PetscFunctionReturn(0);
36250f6de3fSJed Brown }
36350f6de3fSJed Brown 
36404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
36504d965bbSLois Curfman McInnes /*
36604d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
36704d965bbSLois Curfman McInnes 
36804d965bbSLois Curfman McInnes    Input Parameter:
36904d965bbSLois Curfman McInnes .  snes - the SNES context
37004d965bbSLois Curfman McInnes 
37104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
37204d965bbSLois Curfman McInnes */
3734a2ae208SSatish Balay #undef __FUNCT__
37404d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONLS"
37504d7464bSBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes)
3765e42470aSBarry Smith {
377dfbe8321SBarry Smith   PetscErrorCode ierr;
378f1c6b773SPeter Brune   SNESLineSearch linesearch;
3795e42470aSBarry Smith 
3803a40ed3dSBarry Smith   PetscFunctionBegin;
38104d7464bSBarry Smith   ierr = PetscOptionsHead("SNESNEWTONLS options");CHKERRQ(ierr);
382b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3839e764e56SPeter Brune   /* set the default line search type */
3849e764e56SPeter Brune   if (!snes->linesearch) {
385f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
3861a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr);
3879e764e56SPeter Brune   }
3883a40ed3dSBarry Smith   PetscFunctionReturn(0);
3895e42470aSBarry Smith }
3904827ddcaSBarry Smith 
39104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
392ebe3b25bSBarry Smith /*MC
39304d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
39404d965bbSLois Curfman McInnes 
395ebe3b25bSBarry Smith    Options Database:
3965a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
3975a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
3985a9b6599SPeter Brune .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
39982d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
40082d26c24SPeter 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)
4015a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
40282d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
40382d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
404acbee50cSBarry Smith 
405ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
40604d965bbSLois Curfman McInnes 
407ee3001cbSBarry Smith    Level: beginner
408ee3001cbSBarry Smith 
40904d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder()
4105a9b6599SPeter Brune            SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms()
411ebe3b25bSBarry Smith 
412ebe3b25bSBarry Smith M*/
4134a2ae208SSatish Balay #undef __FUNCT__
41404d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONLS"
415*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
4165e42470aSBarry Smith {
417dfbe8321SBarry Smith   PetscErrorCode ierr;
41804d7464bSBarry Smith   SNES_NEWTONLS  *neP;
4195e42470aSBarry Smith 
4203a40ed3dSBarry Smith   PetscFunctionBegin;
42104d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
42204d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
42304d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
42404d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
42504d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
42604d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
4275e42470aSBarry Smith 
42842f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
4296b6c8463SPeter Brune   snes->usespc  = PETSC_TRUE;
43004d7464bSBarry Smith   ierr          = PetscNewLog(snes,SNES_NEWTONLS,&neP);CHKERRQ(ierr);
4315e42470aSBarry Smith   snes->data    = (void*)neP;
4323a40ed3dSBarry Smith   PetscFunctionReturn(0);
4335e42470aSBarry Smith }
434