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 */ 1515e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1525e42470aSBarry Smith G = snes->work[1]; 1535e42470aSBarry Smith W = snes->work[2]; 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); 159*e4ed7901SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 160*e4ed7901SPeter 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 } 167*e4ed7901SPeter Brune } else { 168*e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 169*e4ed7901SPeter Brune } 170*e4ed7901SPeter 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"); 174*e4ed7901SPeter Brune } else { 175*e4ed7901SPeter Brune fnorm = snes->norm_init; 176*e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 177*e4ed7901SPeter 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: 218ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 219e68848bdSBarry Smith and evaluate G = function(Y) (depends on the line search). 220ea4d3ed3SLois Curfman McInnes */ 22185385478SLisandro Dalcin ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 2221936c542SPeter Brune gnorm = fnorm; 223f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 224f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(linesearch, &lssucceed);CHKERRQ(ierr); 225f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 2261936c542SPeter 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); 2274a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 2281936c542SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 2291936c542SPeter Brune if (domainerror) { 2304936397dSBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2314936397dSBarry Smith PetscFunctionReturn(0); 2324936397dSBarry Smith } 233a7cc72afSBarry Smith if (!lssucceed) { 23450ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 235ace3abfcSBarry Smith PetscBool ismin; 236647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2371e2582c4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 2388a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2393505fcc1SBarry Smith break; 2403505fcc1SBarry Smith } 24150ffb88aSMatthew Knepley } 24285385478SLisandro Dalcin /* Monitor convergence */ 24385385478SLisandro Dalcin ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 24485385478SLisandro Dalcin snes->iter = i+1; 24585385478SLisandro Dalcin snes->norm = fnorm; 24685385478SLisandro Dalcin ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 24785385478SLisandro Dalcin SNESLogConvHistory(snes,snes->norm,lits); 2487a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 24985385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 25085385478SLisandro Dalcin if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 251e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2523f149594SLisandro Dalcin if (snes->reason) break; 25316c95cacSBarry Smith } 25452392280SLois Curfman McInnes if (i == maxits) { 255ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 25685385478SLisandro Dalcin if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 25752392280SLois Curfman McInnes } 2583a40ed3dSBarry Smith PetscFunctionReturn(0); 2595e76c431SBarry Smith } 26004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 26104d965bbSLois Curfman McInnes /* 2624b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2634b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 26404d965bbSLois Curfman McInnes 26504d965bbSLois Curfman McInnes Input Parameter: 26604d965bbSLois Curfman McInnes . snes - the SNES context 26704d965bbSLois Curfman McInnes . x - the solution vector 26804d965bbSLois Curfman McInnes 26904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 27004d965bbSLois Curfman McInnes 27104d965bbSLois Curfman McInnes Notes: 2724b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 27304d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 27404d965bbSLois Curfman McInnes the call to SNESSolve(). 27504d965bbSLois Curfman McInnes */ 2764a2ae208SSatish Balay #undef __FUNCT__ 2774b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 278dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2795e76c431SBarry Smith { 280dfbe8321SBarry Smith PetscErrorCode ierr; 2813a40ed3dSBarry Smith 2823a40ed3dSBarry Smith PetscFunctionBegin; 28358c9b817SLisandro Dalcin ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 2846cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 2851936c542SPeter Brune 2863a40ed3dSBarry Smith PetscFunctionReturn(0); 2875e76c431SBarry Smith } 28804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2896b8b9a38SLisandro Dalcin 2906b8b9a38SLisandro Dalcin #undef __FUNCT__ 2916b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS" 2926b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes) 2936b8b9a38SLisandro Dalcin { 2946b8b9a38SLisandro Dalcin PetscFunctionBegin; 2956b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2966b8b9a38SLisandro Dalcin } 2976b8b9a38SLisandro Dalcin 29804d965bbSLois Curfman McInnes /* 2994b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 3004b27c08aSLois Curfman McInnes with SNESCreate_LS(). 30104d965bbSLois Curfman McInnes 30204d965bbSLois Curfman McInnes Input Parameter: 30304d965bbSLois Curfman McInnes . snes - the SNES context 30404d965bbSLois Curfman McInnes 305de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30604d965bbSLois Curfman McInnes */ 3074a2ae208SSatish Balay #undef __FUNCT__ 3084b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 309dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 3105e76c431SBarry Smith { 311dfbe8321SBarry Smith PetscErrorCode ierr; 3123a40ed3dSBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 3146b8b9a38SLisandro Dalcin ierr = SNESReset_LS(snes);CHKERRQ(ierr); 3155c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 3163a40ed3dSBarry Smith PetscFunctionReturn(0); 3175e76c431SBarry Smith } 31804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 31994298653SBarry Smith 320329e7e40SMatthew Knepley /* 3214b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 32204d965bbSLois Curfman McInnes 32304d965bbSLois Curfman McInnes Input Parameters: 32404d965bbSLois Curfman McInnes . SNES - the SNES context 32504d965bbSLois Curfman McInnes . viewer - visualization context 32604d965bbSLois Curfman McInnes 32704d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 32804d965bbSLois Curfman McInnes */ 3294a2ae208SSatish Balay #undef __FUNCT__ 3304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 3316849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 332a935fc98SLois Curfman McInnes { 333dfbe8321SBarry Smith PetscErrorCode ierr; 334ace3abfcSBarry Smith PetscBool iascii; 335a935fc98SLois Curfman McInnes 3363a40ed3dSBarry Smith PetscFunctionBegin; 3372692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 33832077d6dSBarry Smith if (iascii) { 33950f6de3fSJed Brown } 34050f6de3fSJed Brown PetscFunctionReturn(0); 34150f6de3fSJed Brown } 34250f6de3fSJed Brown 34304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 34404d965bbSLois Curfman McInnes /* 3454b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 34604d965bbSLois Curfman McInnes 34704d965bbSLois Curfman McInnes Input Parameter: 34804d965bbSLois Curfman McInnes . snes - the SNES context 34904d965bbSLois Curfman McInnes 35004d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 35104d965bbSLois Curfman McInnes */ 3524a2ae208SSatish Balay #undef __FUNCT__ 3534b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 3546849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 3555e42470aSBarry Smith { 356dfbe8321SBarry Smith PetscErrorCode ierr; 357f1c6b773SPeter Brune SNESLineSearch linesearch; 3585e42470aSBarry Smith 3593a40ed3dSBarry Smith PetscFunctionBegin; 3609e764e56SPeter Brune ierr = PetscOptionsHead("SNESLS options");CHKERRQ(ierr); 361b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3629e764e56SPeter Brune /* set the default line search type */ 3639e764e56SPeter Brune if (!snes->linesearch) { 364f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 365f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_BT);CHKERRQ(ierr); 3669e764e56SPeter Brune } 3673a40ed3dSBarry Smith PetscFunctionReturn(0); 3685e42470aSBarry Smith } 3694827ddcaSBarry Smith 37004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 371ebe3b25bSBarry Smith /*MC 372ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 37304d965bbSLois Curfman McInnes 374ebe3b25bSBarry Smith Options Database: 375ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 37655d4ea13SBarry Smith . -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 377e106eecfSBarry Smith . -snes_ls_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) 378acbee50cSBarry Smith . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 37955d4ea13SBarry Smith . -snes_ls_monitor - print information about progress of line searches 38055d4ea13SBarry Smith - -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms 381acbee50cSBarry Smith 38204d965bbSLois Curfman McInnes 383ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 38404d965bbSLois Curfman McInnes 385ee3001cbSBarry Smith Level: beginner 386ee3001cbSBarry Smith 38717bae607SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 38817bae607SBarry Smith SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 389b3dd4ab5SBarry Smith SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 390ebe3b25bSBarry Smith 391ebe3b25bSBarry Smith M*/ 392fb2e594dSBarry Smith EXTERN_C_BEGIN 3934a2ae208SSatish Balay #undef __FUNCT__ 3944b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 3957087cfbeSBarry Smith PetscErrorCode SNESCreate_LS(SNES snes) 3965e42470aSBarry Smith { 397dfbe8321SBarry Smith PetscErrorCode ierr; 3984b27c08aSLois Curfman McInnes SNES_LS *neP; 3995e42470aSBarry Smith 4003a40ed3dSBarry Smith PetscFunctionBegin; 401e7788613SBarry Smith snes->ops->setup = SNESSetUp_LS; 402e7788613SBarry Smith snes->ops->solve = SNESSolve_LS; 403e7788613SBarry Smith snes->ops->destroy = SNESDestroy_LS; 404e7788613SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_LS; 405e7788613SBarry Smith snes->ops->view = SNESView_LS; 4066b8b9a38SLisandro Dalcin snes->ops->reset = SNESReset_LS; 4075e42470aSBarry Smith 40842f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 40942f4f86dSBarry Smith snes->usespc = PETSC_FALSE; 41038f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 4115e42470aSBarry Smith snes->data = (void*)neP; 4123a40ed3dSBarry Smith PetscFunctionReturn(0); 4135e42470aSBarry Smith } 414fb2e594dSBarry Smith EXTERN_C_END 415