1c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h> 25e42470aSBarry Smith 38a5d9ee4SBarry Smith /* 4420bcc1bSBarry Smith This file implements a truncated Newton method with a line search, 5420bcc1bSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 6420bcc1bSBarry Smith and Mat interfaces for linear solvers, vectors, and matrices, 7420bcc1bSBarry Smith respectively. 8420bcc1bSBarry Smith 9420bcc1bSBarry Smith The following basic routines are required for each nonlinear solver: 10420bcc1bSBarry Smith SNESCreate_XXX() - Creates a nonlinear solver context 11420bcc1bSBarry Smith SNESSetFromOptions_XXX() - Sets runtime options 12420bcc1bSBarry Smith SNESSolve_XXX() - Solves the nonlinear system 13420bcc1bSBarry Smith SNESDestroy_XXX() - Destroys the nonlinear solver context 14420bcc1bSBarry Smith The suffix "_XXX" denotes a particular implementation, in this case 15420bcc1bSBarry Smith we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 16420bcc1bSBarry Smith systems of nonlinear equations with a line search (LS) method. 17420bcc1bSBarry Smith These routines are actually called via the common user interface 18420bcc1bSBarry Smith routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 19420bcc1bSBarry Smith SNESDestroy(), so the application code interface remains identical 20420bcc1bSBarry Smith for all nonlinear solvers. 21420bcc1bSBarry Smith 22420bcc1bSBarry Smith Another key routine is: 23420bcc1bSBarry Smith SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 24420bcc1bSBarry Smith by setting data structures and options. The interface routine SNESSetUp() 25420bcc1bSBarry Smith is not usually called directly by the user, but instead is called by 26420bcc1bSBarry Smith SNESSolve() if necessary. 27420bcc1bSBarry Smith 28420bcc1bSBarry Smith Additional basic routines are: 29420bcc1bSBarry Smith SNESView_XXX() - Prints details of runtime options that 30420bcc1bSBarry Smith have actually been used. 31420bcc1bSBarry Smith These are called by application codes via the interface routines 32420bcc1bSBarry Smith SNESView(). 33420bcc1bSBarry Smith 34420bcc1bSBarry Smith The various types of solvers (preconditioners, Krylov subspace methods, 35420bcc1bSBarry Smith nonlinear solvers, timesteppers) are all organized similarly, so the 36420bcc1bSBarry Smith above description applies to these categories also. 37420bcc1bSBarry Smith 38420bcc1bSBarry Smith */ 39420bcc1bSBarry Smith 40420bcc1bSBarry Smith /* 4120f52e01SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 4220f52e01SBarry 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 4336669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 4420f52e01SBarry Smith for this trick. One assumes that the probability that W is in the null space of J is very, very small. 458a5d9ee4SBarry Smith */ 46d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin) 47d71ae5a4SJacob Faibussowitsch { 488a5d9ee4SBarry Smith PetscReal a1; 49ace3abfcSBarry Smith PetscBool hastranspose; 50f6c5ca78SBarry Smith Vec W; 51728696cfSStefano Zampini PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 528a5d9ee4SBarry Smith 538a5d9ee4SBarry Smith PetscFunctionBegin; 548a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 55728696cfSStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 56728696cfSStefano Zampini if (!objective) { 579566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W)); 5936669109SBarry Smith if (hastranspose) { 608a5d9ee4SBarry Smith /* Compute || J^T F|| */ 619566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, F, W)); 629566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &a1)); 639566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm))); 648a5d9ee4SBarry Smith if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 6536669109SBarry Smith } else { 6636669109SBarry Smith Vec work; 67ea709b57SSatish Balay PetscScalar result; 6836669109SBarry Smith PetscReal wnorm; 6936669109SBarry Smith 709566063dSJacob Faibussowitsch PetscCall(VecSetRandom(W, NULL)); 719566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &work)); 739566063dSJacob Faibussowitsch PetscCall(MatMult(A, W, work)); 749566063dSJacob Faibussowitsch PetscCall(VecDot(F, work, &result)); 759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 7636669109SBarry Smith a1 = PetscAbsScalar(result) / (fnorm * wnorm); 779566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1)); 7836669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 7936669109SBarry Smith } 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W)); 81728696cfSStefano Zampini } 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 838a5d9ee4SBarry Smith } 848a5d9ee4SBarry Smith 8574637425SBarry Smith /* 865ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 8774637425SBarry Smith */ 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X) 89d71ae5a4SJacob Faibussowitsch { 9074637425SBarry Smith PetscReal a1, a2; 91ace3abfcSBarry Smith PetscBool hastranspose; 92b13f2d03SStefano Zampini PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 9374637425SBarry Smith 9474637425SBarry Smith PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 96b13f2d03SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 97b13f2d03SStefano Zampini if (hastranspose && !objective) { 98f6c5ca78SBarry Smith Vec W1, W2; 99f6c5ca78SBarry Smith 1009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W1)); 1019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W2)); 1029566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, W1)); 1039566063dSJacob Faibussowitsch PetscCall(VecAXPY(W1, -1.0, F)); 10474637425SBarry Smith 10574637425SBarry Smith /* Compute || J^T W|| */ 1069566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 1079566063dSJacob Faibussowitsch PetscCall(VecNorm(W1, NORM_2, &a1)); 1089566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &a2)); 10948a46eb9SPierre 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))); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W1)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2)); 11274637425SBarry Smith } 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11474637425SBarry Smith } 11574637425SBarry Smith 116ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars 1175e42470aSBarry Smith /* 11804d7464bSBarry Smith SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 11904d965bbSLois Curfman McInnes method with a line search. 1205e76c431SBarry Smith 1212fe279fdSBarry Smith Input Parameter: 12204d965bbSLois Curfman McInnes . snes - the SNES context 1235e76c431SBarry Smith 1245e76c431SBarry Smith */ 12566976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 126d71ae5a4SJacob Faibussowitsch { 127a7cc72afSBarry Smith PetscInt maxits, i, lits; 128422a814eSBarry Smith SNESLineSearchReason lssucceed; 1294279555eSSatish Balay PetscReal fnorm, xnorm, ynorm; 130f6c5ca78SBarry Smith Vec Y, X, F; 131f1c6b773SPeter Brune SNESLineSearch linesearch; 132c40d0f55SPeter Brune SNESConvergedReason reason; 133*a0254a93SStefano Zampini PC pc; 1344279555eSSatish Balay #if defined(PETSC_USE_INFO) 1354279555eSSatish Balay PetscReal gnorm; 1364279555eSSatish Balay #endif 1375e76c431SBarry Smith 1383a40ed3dSBarry Smith PetscFunctionBegin; 1390b121fc5SBarry 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); 140c579b300SPatrick Farrell 1413d4c4710SBarry Smith snes->numFailures = 0; 1423d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 143184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 144184914b5SBarry Smith 1455e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1465e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14739e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 148d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1495e76c431SBarry Smith 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1514c49b128SBarry Smith snes->iter = 0; 15275567043SBarry Smith snes->norm = 0.0; 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1549566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 155ce5a860aSPeter Brune 156ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 157efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1589566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1599566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 160eee3b80eSStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 1615f3c5e7aSBarry Smith PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER)); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 163b7281c8aSPeter Brune } 164b7281c8aSPeter Brune 1659566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 1669566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 167ce5a860aSPeter Brune } else { 168e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1699566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 1701aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 171c1c75074SPeter Brune } 1721aa26658SKarl Rupp 1739566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 174422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1765e42470aSBarry Smith snes->norm = fnorm; 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1789566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 1793f149594SLisandro Dalcin 18085385478SLisandro Dalcin /* test convergence */ 181d76a863bSStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 1822d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 1833ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 184d034289bSBarry Smith 185*a0254a93SStefano Zampini /* hook state vector to BFGS preconditioner */ 186*a0254a93SStefano Zampini PetscCall(KSPGetPC(snes->ksp, &pc)); 187*a0254a93SStefano Zampini PetscCall(PCLMVMSetUpdateVec(pc, X)); 188*a0254a93SStefano Zampini 1895e76c431SBarry Smith for (i = 0; i < maxits; i++) { 190329e7e40SMatthew Knepley /* Call general purpose update function */ 191dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 192329e7e40SMatthew Knepley 193ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 194efd4aadfSBarry Smith if (snes->npc) { 195efd4aadfSBarry Smith if (snes->npcside == PC_RIGHT) { 1969566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F)); 1979566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1989566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 1999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 2009566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 201eee3b80eSStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 2025f3c5e7aSBarry Smith PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER)); 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 204c40d0f55SPeter Brune } 2059566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 206efd4aadfSBarry Smith } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2079566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, F)); 2089566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 209eee3b80eSStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 2105f3c5e7aSBarry Smith PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER)); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212b7281c8aSPeter Brune } 213ce5a860aSPeter Brune } 214c40d0f55SPeter Brune } 215c40d0f55SPeter Brune 216ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 2179566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 21807b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2199566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2209566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 221422a814eSBarry Smith SNESCheckKSPSolve(snes); 2229566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 22363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 22474637425SBarry Smith 2251baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y)); 22674637425SBarry Smith 2274279555eSSatish Balay #if defined(PETSC_USE_INFO) 2284279555eSSatish Balay gnorm = fnorm; 2294279555eSSatish Balay #endif 230ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 231d12b9ff3SPeter Brune X <- X - lambda*Y 232d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 233ea4d3ed3SLois Curfman McInnes */ 2349566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 2359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed)); 2369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 2379566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed)); 2387e49c775SBarry Smith if (snes->reason) break; 239422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 240422a814eSBarry Smith if (lssucceed) { 241c21ba15cSPeter Brune if (snes->stol * xnorm > ynorm) { 242c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 244c21ba15cSPeter Brune } 24550ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 246ace3abfcSBarry Smith PetscBool ismin; 247647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2489566063dSJacob Faibussowitsch PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin)); 2498a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2503add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2513add74b1SBarry Smith PetscViewer monitor; 2529566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 25328b400f6SJacob 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]); 254f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]); 2553add74b1SBarry Smith } 2563505fcc1SBarry Smith break; 2573505fcc1SBarry Smith } 25850ffb88aSMatthew Knepley } 25985385478SLisandro Dalcin /* Monitor convergence */ 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 26185385478SLisandro Dalcin snes->iter = i + 1; 26285385478SLisandro Dalcin snes->norm = fnorm; 263c1e67a49SFande Kong snes->ynorm = ynorm; 264c1e67a49SFande Kong snes->xnorm = xnorm; 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2669566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 2672a997342SPeter Brune /* Test for convergence */ 2682d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 2692d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2703f149594SLisandro Dalcin if (snes->reason) break; 27116c95cacSBarry Smith } 2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2735e76c431SBarry Smith } 274f6dfbefdSBarry Smith 27504d965bbSLois Curfman McInnes /* 27604d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 27704d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 27804d965bbSLois Curfman McInnes 27904d965bbSLois Curfman McInnes Input Parameter: 28004d965bbSLois Curfman McInnes . snes - the SNES context 28104d965bbSLois Curfman McInnes . x - the solution vector 28204d965bbSLois Curfman McInnes 28304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 28404d965bbSLois Curfman McInnes 28504d965bbSLois Curfman McInnes */ 28666976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 287d71ae5a4SJacob Faibussowitsch { 2883a40ed3dSBarry Smith PetscFunctionBegin; 2899566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 290efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 2913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2925e76c431SBarry Smith } 2936b8b9a38SLisandro Dalcin 29466976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NEWTONLS(SNES snes) 295d71ae5a4SJacob Faibussowitsch { 2966b8b9a38SLisandro Dalcin PetscFunctionBegin; 2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2986b8b9a38SLisandro Dalcin } 2996b8b9a38SLisandro Dalcin 30004d965bbSLois Curfman McInnes /* 30104d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 30204d7464bSBarry Smith with SNESCreate_NEWTONLS(). 30304d965bbSLois Curfman McInnes 30404d965bbSLois Curfman McInnes Input Parameter: 30504d965bbSLois Curfman McInnes . snes - the SNES context 30604d965bbSLois Curfman McInnes 307de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30804d965bbSLois Curfman McInnes */ 30966976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 310d71ae5a4SJacob Faibussowitsch { 3113a40ed3dSBarry Smith PetscFunctionBegin; 3129566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONLS(snes)); 3139566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3155e76c431SBarry Smith } 31694298653SBarry Smith 317329e7e40SMatthew Knepley /* 31804d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 31904d965bbSLois Curfman McInnes 32004d965bbSLois Curfman McInnes Input Parameters: 32104d965bbSLois Curfman McInnes . SNES - the SNES context 32204d965bbSLois Curfman McInnes . viewer - visualization context 32304d965bbSLois Curfman McInnes 32404d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 32504d965bbSLois Curfman McInnes */ 326d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) 327d71ae5a4SJacob Faibussowitsch { 328ace3abfcSBarry Smith PetscBool iascii; 329a935fc98SLois Curfman McInnes 3303a40ed3dSBarry Smith PetscFunctionBegin; 3319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3329371c9d4SSatish Balay if (iascii) { } 3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33450f6de3fSJed Brown } 33550f6de3fSJed Brown 33604d965bbSLois Curfman McInnes /* 33704d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 33804d965bbSLois Curfman McInnes 33904d965bbSLois Curfman McInnes Input Parameter: 34004d965bbSLois Curfman McInnes . snes - the SNES context 34104d965bbSLois Curfman McInnes 34204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 34304d965bbSLois Curfman McInnes */ 344d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) 345d71ae5a4SJacob Faibussowitsch { 3463a40ed3dSBarry Smith PetscFunctionBegin; 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3485e42470aSBarry Smith } 3494827ddcaSBarry Smith 350ebe3b25bSBarry Smith /*MC 35104d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 35204d965bbSLois Curfman McInnes 3533c7db156SBarry Smith Options Database Keys: 3545a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3555a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 356f6dfbefdSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`) 35782d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 35882d26c24SPeter 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) 3595a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 36082d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 36182d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 362acbee50cSBarry Smith 363420bcc1bSBarry Smith Level: beginner 364420bcc1bSBarry Smith 365f6dfbefdSBarry Smith Note: 366f6dfbefdSBarry Smith This is the default nonlinear solver in `SNES` 36704d965bbSLois Curfman McInnes 368420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()` 369f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()` 370ebe3b25bSBarry Smith M*/ 371d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 372d71ae5a4SJacob Faibussowitsch { 37304d7464bSBarry Smith SNES_NEWTONLS *neP; 374d8d34be6SBarry Smith SNESLineSearch linesearch; 3755e42470aSBarry Smith 3763a40ed3dSBarry Smith PetscFunctionBegin; 37704d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 37804d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 37904d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 38004d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 38104d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 38204d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 3835e42470aSBarry Smith 384efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 38542f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 386efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 3874fc747eaSLawrence Mitchell 3889566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 38948a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 390d8d34be6SBarry Smith 3914fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3924fc747eaSLawrence Mitchell 3934dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 3945e42470aSBarry Smith snes->data = (void *)neP; 3953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3965e42470aSBarry Smith } 397