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__ 114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 12ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_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 32abb0e124SMatthew Knepley ierr = VecSetRandom(W,PETSC_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__ 494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 501e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_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 864b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) 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 /* 1114b27c08aSLois Curfman McInnes SNESSolve_LS - 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__ 1304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 131dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(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; 13885385478SLisandro Dalcin Vec Y,X,F,G,W; 1393d4c4710SBarry Smith KSPConvergedReason kspreason; 1401936c542SPeter Brune PetscBool domainerror; 141f1c6b773SPeter Brune SNESLineSearch linesearch; 1425e76c431SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 1443d4c4710SBarry Smith snes->numFailures = 0; 1453d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 146184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 147184914b5SBarry Smith 1485e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1495e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 15039e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 151d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 152d12b9ff3SPeter Brune G = snes->work[0]; 153d12b9ff3SPeter Brune W = snes->work[1]; 1545e76c431SBarry Smith 1554c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1564c49b128SBarry Smith snes->iter = 0; 15775567043SBarry Smith snes->norm = 0.0; 1584c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 159e4ed7901SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 160e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 16119717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1621936c542SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1631936c542SPeter Brune if (domainerror) { 1644936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1654936397dSBarry Smith PetscFunctionReturn(0); 1664936397dSBarry Smith } 167e4ed7901SPeter Brune } else { 168e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 169e4ed7901SPeter Brune } 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); 17362cbcd01SMatthew G Knepley if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 174e4ed7901SPeter Brune } else { 175e4ed7901SPeter Brune fnorm = snes->norm_init; 176e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 177e4ed7901SPeter Brune } 1780f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1795e42470aSBarry Smith snes->norm = fnorm; 1800f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 18184c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 1827a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1833f149594SLisandro Dalcin 1843f149594SLisandro Dalcin /* set parameter for default relative tolerance convergence test */ 1853f149594SLisandro Dalcin snes->ttol = fnorm*snes->rtol; 18685385478SLisandro Dalcin /* test convergence */ 18785385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 18806ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 189d034289bSBarry Smith 1905e76c431SBarry Smith for (i=0; i<maxits; i++) { 1915e76c431SBarry Smith 192329e7e40SMatthew Knepley /* Call general purpose update function */ 193e7788613SBarry Smith if (snes->ops->update) { 194e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 195de8cb200SMatthew Knepley } 196329e7e40SMatthew Knepley 197ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1985334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 19994b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 20071f87433Sdalcinl ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 2013d4c4710SBarry Smith ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 2023d4c4710SBarry Smith if (kspreason < 0) { 2033d4c4710SBarry Smith if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 204ef998cc9SBarry 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); 2053d4c4710SBarry Smith snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 206ab81109fSSatish Balay break; 2073d4c4710SBarry Smith } 2083d4c4710SBarry Smith } 2093d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2103f149594SLisandro Dalcin snes->linear_its += lits; 2113f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 21274637425SBarry Smith 213b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 2141e2582c4SBarry Smith ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 21585471664SBarry Smith } 21674637425SBarry Smith 217ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 218d12b9ff3SPeter Brune X <- X - lambda*Y 219d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 220ea4d3ed3SLois Curfman McInnes */ 2211936c542SPeter Brune gnorm = fnorm; 222f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 223f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(linesearch, &lssucceed);CHKERRQ(ierr); 224f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 2251936c542SPeter 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); 2264a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2271936c542SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 2281936c542SPeter Brune if (domainerror) { 2294936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2304936397dSBarry Smith PetscFunctionReturn(0); 2314936397dSBarry Smith } 232a7cc72afSBarry Smith if (!lssucceed) { 23350ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 234ace3abfcSBarry Smith PetscBool ismin; 235647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 236a7a23023SPeter Brune ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,F,X,fnorm,&ismin);CHKERRQ(ierr); 2378a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2383505fcc1SBarry Smith break; 2393505fcc1SBarry Smith } 24050ffb88aSMatthew Knepley } 24185385478SLisandro Dalcin /* Monitor convergence */ 24285385478SLisandro Dalcin ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 24385385478SLisandro Dalcin snes->iter = i+1; 24485385478SLisandro Dalcin snes->norm = fnorm; 24585385478SLisandro Dalcin ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 24685385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 2477a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 24885385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 24985385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 250e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2513f149594SLisandro Dalcin if (snes->reason) break; 25216c95cacSBarry Smith } 25352392280SLois Curfman McInnes if (i == maxits) { 254ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 25585385478SLisandro Dalcin if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 25652392280SLois Curfman McInnes } 2573a40ed3dSBarry Smith PetscFunctionReturn(0); 2585e76c431SBarry Smith } 25904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 26004d965bbSLois Curfman McInnes /* 2614b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2624b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 26304d965bbSLois Curfman McInnes 26404d965bbSLois Curfman McInnes Input Parameter: 26504d965bbSLois Curfman McInnes . snes - the SNES context 26604d965bbSLois Curfman McInnes . x - the solution vector 26704d965bbSLois Curfman McInnes 26804d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 26904d965bbSLois Curfman McInnes 27004d965bbSLois Curfman McInnes Notes: 2714b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 27204d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 27304d965bbSLois Curfman McInnes the call to SNESSolve(). 27404d965bbSLois Curfman McInnes */ 2754a2ae208SSatish Balay #undef __FUNCT__ 2764b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 277dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2785e76c431SBarry Smith { 279dfbe8321SBarry Smith PetscErrorCode ierr; 2803a40ed3dSBarry Smith 2813a40ed3dSBarry Smith PetscFunctionBegin; 282d12b9ff3SPeter Brune ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr); 2836cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 2841936c542SPeter Brune 2853a40ed3dSBarry Smith PetscFunctionReturn(0); 2865e76c431SBarry Smith } 28704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2886b8b9a38SLisandro Dalcin 2896b8b9a38SLisandro Dalcin #undef __FUNCT__ 2906b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS" 2916b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes) 2926b8b9a38SLisandro Dalcin { 2936b8b9a38SLisandro Dalcin PetscFunctionBegin; 2946b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2956b8b9a38SLisandro Dalcin } 2966b8b9a38SLisandro Dalcin 29704d965bbSLois Curfman McInnes /* 2984b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2994b27c08aSLois Curfman McInnes with SNESCreate_LS(). 30004d965bbSLois Curfman McInnes 30104d965bbSLois Curfman McInnes Input Parameter: 30204d965bbSLois Curfman McInnes . snes - the SNES context 30304d965bbSLois Curfman McInnes 304de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30504d965bbSLois Curfman McInnes */ 3064a2ae208SSatish Balay #undef __FUNCT__ 3074b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 308dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3095e76c431SBarry Smith { 310dfbe8321SBarry Smith PetscErrorCode ierr; 3113a40ed3dSBarry Smith 3123a40ed3dSBarry Smith PetscFunctionBegin; 3136b8b9a38SLisandro Dalcin ierr = SNESReset_LS(snes);CHKERRQ(ierr); 3145c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 3165e76c431SBarry Smith } 31704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 31894298653SBarry Smith 319329e7e40SMatthew Knepley /* 3204b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 32104d965bbSLois Curfman McInnes 32204d965bbSLois Curfman McInnes Input Parameters: 32304d965bbSLois Curfman McInnes . SNES - the SNES context 32404d965bbSLois Curfman McInnes . viewer - visualization context 32504d965bbSLois Curfman McInnes 32604d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 32704d965bbSLois Curfman McInnes */ 3284a2ae208SSatish Balay #undef __FUNCT__ 3294b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 3306849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 331a935fc98SLois Curfman McInnes { 332dfbe8321SBarry Smith PetscErrorCode ierr; 333ace3abfcSBarry Smith PetscBool iascii; 334a935fc98SLois Curfman McInnes 3353a40ed3dSBarry Smith PetscFunctionBegin; 336251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 33732077d6dSBarry Smith if (iascii) { 33850f6de3fSJed Brown } 33950f6de3fSJed Brown PetscFunctionReturn(0); 34050f6de3fSJed Brown } 34150f6de3fSJed Brown 34204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 34304d965bbSLois Curfman McInnes /* 3444b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 34504d965bbSLois Curfman McInnes 34604d965bbSLois Curfman McInnes Input Parameter: 34704d965bbSLois Curfman McInnes . snes - the SNES context 34804d965bbSLois Curfman McInnes 34904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 35004d965bbSLois Curfman McInnes */ 3514a2ae208SSatish Balay #undef __FUNCT__ 3524b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 3536849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 3545e42470aSBarry Smith { 355dfbe8321SBarry Smith PetscErrorCode ierr; 356f1c6b773SPeter Brune SNESLineSearch linesearch; 3575e42470aSBarry Smith 3583a40ed3dSBarry Smith PetscFunctionBegin; 3599e764e56SPeter Brune ierr = PetscOptionsHead("SNESLS options");CHKERRQ(ierr); 360b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3619e764e56SPeter Brune /* set the default line search type */ 3629e764e56SPeter Brune if (!snes->linesearch) { 363f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 3641a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 3659e764e56SPeter Brune } 3663a40ed3dSBarry Smith PetscFunctionReturn(0); 3675e42470aSBarry Smith } 3684827ddcaSBarry Smith 36904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 370ebe3b25bSBarry Smith /*MC 371ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 37204d965bbSLois Curfman McInnes 373ebe3b25bSBarry Smith Options Database: 374*5a9b6599SPeter Brune + -snes_linesearch_type<bt> - bt,basic. Select line search type 375*5a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 376*5a9b6599SPeter Brune . -snes_linesearch_norms<true> - Turns on/off computation of the norms for basic linesearch 37782d26c24SPeter Brune . -snes_linesearch_alpha<alpha> - Sets alpha used in determining if reduction in function norm is sufficient 37882d26c24SPeter 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) 379*5a9b6599SPeter Brune . -snes_linesearch_minlambda<minlambda> - Sets the minimum lambda the line search will tolerate 38082d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 38182d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 382acbee50cSBarry Smith 383ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 38404d965bbSLois Curfman McInnes 385ee3001cbSBarry Smith Level: beginner 386ee3001cbSBarry Smith 387*5a9b6599SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder() 388*5a9b6599SPeter Brune SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms() 389ebe3b25bSBarry Smith 390ebe3b25bSBarry Smith M*/ 391fb2e594dSBarry Smith EXTERN_C_BEGIN 3924a2ae208SSatish Balay #undef __FUNCT__ 3934b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 3947087cfbeSBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 3955e42470aSBarry Smith { 396dfbe8321SBarry Smith PetscErrorCode ierr; 3974b27c08aSLois Curfman McInnes SNES_LS *neP; 3985e42470aSBarry Smith 3993a40ed3dSBarry Smith PetscFunctionBegin; 400e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 401e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 402e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 403e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 404e7788613SBarry Smith snes->ops->view = SNESView_LS; 4056b8b9a38SLisandro Dalcin snes->ops->reset = SNESReset_LS; 4065e42470aSBarry Smith 40742f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 40842f4f86dSBarry Smith snes->usespc = PETSC_FALSE; 40938f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 4105e42470aSBarry Smith snes->data = (void*)neP; 4113a40ed3dSBarry Smith PetscFunctionReturn(0); 4125e42470aSBarry Smith } 413fb2e594dSBarry Smith EXTERN_C_END 414