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 */ 109371c9d4SSatish Balay static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin) { 118a5d9ee4SBarry Smith PetscReal a1; 12ace3abfcSBarry Smith PetscBool hastranspose; 13f6c5ca78SBarry Smith Vec W; 148a5d9ee4SBarry Smith 158a5d9ee4SBarry Smith PetscFunctionBegin; 168a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 179566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W)); 1936669109SBarry Smith if (hastranspose) { 208a5d9ee4SBarry Smith /* Compute || J^T F|| */ 219566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, F, W)); 229566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &a1)); 239566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm))); 248a5d9ee4SBarry Smith if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 2536669109SBarry Smith } else { 2636669109SBarry Smith Vec work; 27ea709b57SSatish Balay PetscScalar result; 2836669109SBarry Smith PetscReal wnorm; 2936669109SBarry Smith 309566063dSJacob Faibussowitsch PetscCall(VecSetRandom(W, NULL)); 319566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &work)); 339566063dSJacob Faibussowitsch PetscCall(MatMult(A, W, work)); 349566063dSJacob Faibussowitsch PetscCall(VecDot(F, work, &result)); 359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 3636669109SBarry Smith a1 = PetscAbsScalar(result) / (fnorm * wnorm); 379566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1)); 3836669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 3936669109SBarry Smith } 409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W)); 418a5d9ee4SBarry Smith PetscFunctionReturn(0); 428a5d9ee4SBarry Smith } 438a5d9ee4SBarry Smith 4474637425SBarry Smith /* 455ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4674637425SBarry Smith */ 479371c9d4SSatish Balay static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X) { 4874637425SBarry Smith PetscReal a1, a2; 49ace3abfcSBarry Smith PetscBool hastranspose; 5074637425SBarry Smith 5174637425SBarry Smith PetscFunctionBegin; 529566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 5374637425SBarry Smith if (hastranspose) { 54f6c5ca78SBarry Smith Vec W1, W2; 55f6c5ca78SBarry Smith 569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W1)); 579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W2)); 589566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, W1)); 599566063dSJacob Faibussowitsch PetscCall(VecAXPY(W1, -1.0, F)); 6074637425SBarry Smith 6174637425SBarry Smith /* Compute || J^T W|| */ 629566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 639566063dSJacob Faibussowitsch PetscCall(VecNorm(W1, NORM_2, &a1)); 649566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &a2)); 65*48a46eb9SPierre Jolivet if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n", (double)(a2 / a1))); 669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W1)); 679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2)); 6874637425SBarry Smith } 6974637425SBarry Smith PetscFunctionReturn(0); 7074637425SBarry Smith } 7174637425SBarry Smith 7204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7304d965bbSLois Curfman McInnes 7404d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7594b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7604d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 7704d965bbSLois Curfman McInnes respectively. 7804d965bbSLois Curfman McInnes 79fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8004d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8104d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 82fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8304d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8404d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 8504d7464bSBarry Smith we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 864b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 8704d965bbSLois Curfman McInnes These routines are actually called via the common user interface 8804d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 8904d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9004d965bbSLois Curfman McInnes for all nonlinear solvers. 9104d965bbSLois Curfman McInnes 9204d965bbSLois Curfman McInnes Another key routine is: 9304d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9404d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9504d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9604d965bbSLois Curfman McInnes SNESSolve() if necessary. 9704d965bbSLois Curfman McInnes 9804d965bbSLois Curfman McInnes Additional basic routines are: 9904d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10004d965bbSLois Curfman McInnes have actually been used. 10104d965bbSLois Curfman McInnes These are called by application codes via the interface routines 102186905e3SBarry Smith SNESView(). 10304d965bbSLois Curfman McInnes 10404d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10504d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10604d965bbSLois Curfman McInnes above description applies to these categories also. 10704d965bbSLois Curfman McInnes 10804d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1095e42470aSBarry Smith /* 11004d7464bSBarry Smith SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 11104d965bbSLois Curfman McInnes method with a line search. 1125e76c431SBarry Smith 11304d965bbSLois Curfman McInnes Input Parameters: 11404d965bbSLois Curfman McInnes . snes - the SNES context 1155e76c431SBarry Smith 11604d965bbSLois Curfman McInnes Output Parameter: 117c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 11804d965bbSLois Curfman McInnes 11904d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1205e76c431SBarry Smith 1215e76c431SBarry Smith Notes: 1225e76c431SBarry Smith This implements essentially a truncated Newton method with a 1235e76c431SBarry Smith line search. By default a cubic backtracking line search 1245e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1255e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 126393d2d9aSLois Curfman McInnes and Schnabel. 1275e76c431SBarry Smith */ 1289371c9d4SSatish Balay PetscErrorCode SNESSolve_NEWTONLS(SNES snes) { 129a7cc72afSBarry Smith PetscInt maxits, i, lits; 130422a814eSBarry Smith SNESLineSearchReason lssucceed; 1311936c542SPeter Brune PetscReal fnorm, gnorm, xnorm, ynorm; 132f6c5ca78SBarry Smith Vec Y, X, F; 133f1c6b773SPeter Brune SNESLineSearch linesearch; 134c40d0f55SPeter Brune SNESConvergedReason reason; 1355e76c431SBarry Smith 1363a40ed3dSBarry Smith PetscFunctionBegin; 1370b121fc5SBarry Smith PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 138c579b300SPatrick Farrell 1393d4c4710SBarry Smith snes->numFailures = 0; 1403d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 141184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 142184914b5SBarry Smith 1435e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1445e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14539e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 146d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1475e76c431SBarry Smith 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1494c49b128SBarry Smith snes->iter = 0; 15075567043SBarry Smith snes->norm = 0.0; 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1529566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 153ce5a860aSPeter Brune 154ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 155efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1569566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1579566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 158b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 159b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 160b7281c8aSPeter Brune PetscFunctionReturn(0); 161b7281c8aSPeter Brune } 162b7281c8aSPeter Brune 1639566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 1649566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 165ce5a860aSPeter Brune } else { 166e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1679566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 1681aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 169c1c75074SPeter Brune } 1701aa26658SKarl Rupp 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 172422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1745e42470aSBarry Smith snes->norm = fnorm; 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1769566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 1779566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 1783f149594SLisandro Dalcin 17985385478SLisandro Dalcin /* test convergence */ 180dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 18106ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 182d034289bSBarry Smith 1835e76c431SBarry Smith for (i = 0; i < maxits; i++) { 184329e7e40SMatthew Knepley /* Call general purpose update function */ 185dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 186329e7e40SMatthew Knepley 187ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 188efd4aadfSBarry Smith if (snes->npc) { 189efd4aadfSBarry Smith if (snes->npcside == PC_RIGHT) { 1909566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F)); 1919566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1929566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 1939566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1949566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 195c40d0f55SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 196c40d0f55SPeter Brune snes->reason = SNES_DIVERGED_INNER; 197c40d0f55SPeter Brune PetscFunctionReturn(0); 198c40d0f55SPeter Brune } 1999566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 200efd4aadfSBarry Smith } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2019566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, F)); 2029566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 203b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 204b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 205b7281c8aSPeter Brune PetscFunctionReturn(0); 206b7281c8aSPeter Brune } 207ce5a860aSPeter Brune } 208c40d0f55SPeter Brune } 209c40d0f55SPeter Brune 210ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 2119566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 21207b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2139566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2149566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 215422a814eSBarry Smith SNESCheckKSPSolve(snes); 2169566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 21763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 21874637425SBarry Smith 2191baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y)); 22074637425SBarry Smith 221ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 222d12b9ff3SPeter Brune X <- X - lambda*Y 223d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 224ea4d3ed3SLois Curfman McInnes */ 2251936c542SPeter Brune gnorm = fnorm; 2269566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 2279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed)); 2289566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 2299566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed)); 2307e49c775SBarry Smith if (snes->reason) break; 231422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 232422a814eSBarry Smith if (lssucceed) { 233c21ba15cSPeter Brune if (snes->stol * xnorm > ynorm) { 234c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 235c21ba15cSPeter Brune PetscFunctionReturn(0); 236c21ba15cSPeter Brune } 23750ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 238ace3abfcSBarry Smith PetscBool ismin; 239647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2409566063dSJacob Faibussowitsch PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin)); 2418a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2423add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2433add74b1SBarry Smith PetscViewer monitor; 2449566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 24528b400f6SJacob Faibussowitsch PetscCheck(monitor, PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due to %s. Suggest running with -snes_linesearch_monitor", SNESConvergedReasons[snes->reason]); 246f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]); 2473add74b1SBarry Smith } 2483505fcc1SBarry Smith break; 2493505fcc1SBarry Smith } 25050ffb88aSMatthew Knepley } 25185385478SLisandro Dalcin /* Monitor convergence */ 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 25385385478SLisandro Dalcin snes->iter = i + 1; 25485385478SLisandro Dalcin snes->norm = fnorm; 255c1e67a49SFande Kong snes->ynorm = ynorm; 256c1e67a49SFande Kong snes->xnorm = xnorm; 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2589566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 2599566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2602a997342SPeter Brune /* Test for convergence */ 261dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 2623f149594SLisandro Dalcin if (snes->reason) break; 26316c95cacSBarry Smith } 26452392280SLois Curfman McInnes if (i == maxits) { 26563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 26685385478SLisandro Dalcin if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 26752392280SLois Curfman McInnes } 2683a40ed3dSBarry Smith PetscFunctionReturn(0); 2695e76c431SBarry Smith } 27004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 27104d965bbSLois Curfman McInnes /* 27204d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 27304d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 27404d965bbSLois Curfman McInnes 27504d965bbSLois Curfman McInnes Input Parameter: 27604d965bbSLois Curfman McInnes . snes - the SNES context 27704d965bbSLois Curfman McInnes . x - the solution vector 27804d965bbSLois Curfman McInnes 27904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 28004d965bbSLois Curfman McInnes 28104d965bbSLois Curfman McInnes Notes: 2824b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 28304d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 28404d965bbSLois Curfman McInnes the call to SNESSolve(). 28504d965bbSLois Curfman McInnes */ 2869371c9d4SSatish Balay PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) { 2873a40ed3dSBarry Smith PetscFunctionBegin; 2889566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 289efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 2903a40ed3dSBarry Smith PetscFunctionReturn(0); 2915e76c431SBarry Smith } 29204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2936b8b9a38SLisandro Dalcin 2949371c9d4SSatish Balay PetscErrorCode SNESReset_NEWTONLS(SNES snes) { 2956b8b9a38SLisandro Dalcin PetscFunctionBegin; 2966b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2976b8b9a38SLisandro Dalcin } 2986b8b9a38SLisandro Dalcin 29904d965bbSLois Curfman McInnes /* 30004d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 30104d7464bSBarry Smith with SNESCreate_NEWTONLS(). 30204d965bbSLois Curfman McInnes 30304d965bbSLois Curfman McInnes Input Parameter: 30404d965bbSLois Curfman McInnes . snes - the SNES context 30504d965bbSLois Curfman McInnes 306de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30704d965bbSLois Curfman McInnes */ 3089371c9d4SSatish Balay PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) { 3093a40ed3dSBarry Smith PetscFunctionBegin; 3109566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONLS(snes)); 3119566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3123a40ed3dSBarry Smith PetscFunctionReturn(0); 3135e76c431SBarry Smith } 31404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 31594298653SBarry Smith 316329e7e40SMatthew Knepley /* 31704d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 31804d965bbSLois Curfman McInnes 31904d965bbSLois Curfman McInnes Input Parameters: 32004d965bbSLois Curfman McInnes . SNES - the SNES context 32104d965bbSLois Curfman McInnes . viewer - visualization context 32204d965bbSLois Curfman McInnes 32304d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 32404d965bbSLois Curfman McInnes */ 3259371c9d4SSatish Balay static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) { 326ace3abfcSBarry Smith PetscBool iascii; 327a935fc98SLois Curfman McInnes 3283a40ed3dSBarry Smith PetscFunctionBegin; 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3309371c9d4SSatish Balay if (iascii) { } 33150f6de3fSJed Brown PetscFunctionReturn(0); 33250f6de3fSJed Brown } 33350f6de3fSJed Brown 33404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 33504d965bbSLois Curfman McInnes /* 33604d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 33704d965bbSLois Curfman McInnes 33804d965bbSLois Curfman McInnes Input Parameter: 33904d965bbSLois Curfman McInnes . snes - the SNES context 34004d965bbSLois Curfman McInnes 34104d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 34204d965bbSLois Curfman McInnes */ 3439371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) { 3443a40ed3dSBarry Smith PetscFunctionBegin; 3453a40ed3dSBarry Smith PetscFunctionReturn(0); 3465e42470aSBarry Smith } 3474827ddcaSBarry Smith 34804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 349ebe3b25bSBarry Smith /*MC 35004d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 35104d965bbSLois Curfman McInnes 352ebe3b25bSBarry Smith Options Database: 3535a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3545a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 3551fe24845SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (SNESLineSearchSetComputeNorms()) 35682d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 35782d26c24SPeter 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) 3585a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 35982d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 36082d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 361acbee50cSBarry Smith 36295452b02SPatrick Sanan Notes: 36395452b02SPatrick Sanan This is the default nonlinear solver in SNES 36404d965bbSLois Curfman McInnes 365ee3001cbSBarry Smith Level: beginner 366ee3001cbSBarry Smith 367db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()` 368db781477SPatrick Sanan `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()` 369ebe3b25bSBarry Smith 370ebe3b25bSBarry Smith M*/ 3719371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) { 37204d7464bSBarry Smith SNES_NEWTONLS *neP; 373d8d34be6SBarry Smith SNESLineSearch linesearch; 3745e42470aSBarry Smith 3753a40ed3dSBarry Smith PetscFunctionBegin; 37604d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 37704d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 37804d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 37904d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 38004d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 38104d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 3825e42470aSBarry Smith 383efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 38442f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 385efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 3864fc747eaSLawrence Mitchell 3879566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 388*48a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 389d8d34be6SBarry Smith 3904fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3914fc747eaSLawrence Mitchell 3929566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes, &neP)); 3935e42470aSBarry Smith snes->data = (void *)neP; 3943a40ed3dSBarry Smith PetscFunctionReturn(0); 3955e42470aSBarry Smith } 396