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; 133a0254a93SStefano 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 185a0254a93SStefano Zampini /* hook state vector to BFGS preconditioner */ 186a0254a93SStefano Zampini PetscCall(KSPGetPC(snes->ksp, &pc)); 187a0254a93SStefano Zampini PetscCall(PCLMVMSetUpdateVec(pc, X)); 188a0254a93SStefano Zampini 1895e76c431SBarry Smith for (i = 0; i < maxits; i++) { 190329e7e40SMatthew Knepley /* Call general purpose update function */ 191dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 192*3238a83eSBarry Smith PetscCall(VecNorm(snes->vec_func, NORM_2, &fnorm)); /* no-op unless update() function changed f() */ 193329e7e40SMatthew Knepley 194ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 195efd4aadfSBarry Smith if (snes->npc) { 196efd4aadfSBarry Smith if (snes->npcside == PC_RIGHT) { 1979566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F)); 1989566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1999566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 2009566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 2019566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 202eee3b80eSStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 2035f3c5e7aSBarry Smith PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER)); 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 205c40d0f55SPeter Brune } 2069566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 207efd4aadfSBarry Smith } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2089566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, F)); 2099566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 210eee3b80eSStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 2115f3c5e7aSBarry Smith PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER)); 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 213b7281c8aSPeter Brune } 214ce5a860aSPeter Brune } 215c40d0f55SPeter Brune } 216c40d0f55SPeter Brune 217ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 2189566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 21907b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2209566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2219566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 222422a814eSBarry Smith SNESCheckKSPSolve(snes); 2239566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 22463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 22574637425SBarry Smith 2261baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y)); 22774637425SBarry Smith 2284279555eSSatish Balay #if defined(PETSC_USE_INFO) 2294279555eSSatish Balay gnorm = fnorm; 2304279555eSSatish Balay #endif 231ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 232d12b9ff3SPeter Brune X <- X - lambda*Y 233d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 234ea4d3ed3SLois Curfman McInnes */ 2359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 2369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed)); 2379566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 2389566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed)); 2397e49c775SBarry Smith if (snes->reason) break; 240422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 241422a814eSBarry Smith if (lssucceed) { 242c21ba15cSPeter Brune if (snes->stol * xnorm > ynorm) { 243c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245c21ba15cSPeter Brune } 24650ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 247ace3abfcSBarry Smith PetscBool ismin; 248647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2499566063dSJacob Faibussowitsch PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin)); 2508a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2513add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2523add74b1SBarry Smith PetscViewer monitor; 2539566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 25428b400f6SJacob 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]); 255f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]); 2563add74b1SBarry Smith } 2573505fcc1SBarry Smith break; 2583505fcc1SBarry Smith } 25950ffb88aSMatthew Knepley } 26085385478SLisandro Dalcin /* Monitor convergence */ 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 26285385478SLisandro Dalcin snes->iter = i + 1; 26385385478SLisandro Dalcin snes->norm = fnorm; 264c1e67a49SFande Kong snes->ynorm = ynorm; 265c1e67a49SFande Kong snes->xnorm = xnorm; 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2679566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 2682a997342SPeter Brune /* Test for convergence */ 2692d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 2702d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2713f149594SLisandro Dalcin if (snes->reason) break; 27216c95cacSBarry Smith } 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2745e76c431SBarry Smith } 275f6dfbefdSBarry Smith 27604d965bbSLois Curfman McInnes /* 27704d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 27804d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 27904d965bbSLois Curfman McInnes 28004d965bbSLois Curfman McInnes Input Parameter: 28104d965bbSLois Curfman McInnes . snes - the SNES context 28204d965bbSLois Curfman McInnes . x - the solution vector 28304d965bbSLois Curfman McInnes 28404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 28504d965bbSLois Curfman McInnes 28604d965bbSLois Curfman McInnes */ 28766976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 288d71ae5a4SJacob Faibussowitsch { 2893a40ed3dSBarry Smith PetscFunctionBegin; 2909566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 291efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2935e76c431SBarry Smith } 2946b8b9a38SLisandro Dalcin 29566976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NEWTONLS(SNES snes) 296d71ae5a4SJacob Faibussowitsch { 2976b8b9a38SLisandro Dalcin PetscFunctionBegin; 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2996b8b9a38SLisandro Dalcin } 3006b8b9a38SLisandro Dalcin 30104d965bbSLois Curfman McInnes /* 30204d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 30304d7464bSBarry Smith with SNESCreate_NEWTONLS(). 30404d965bbSLois Curfman McInnes 30504d965bbSLois Curfman McInnes Input Parameter: 30604d965bbSLois Curfman McInnes . snes - the SNES context 30704d965bbSLois Curfman McInnes 308de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30904d965bbSLois Curfman McInnes */ 31066976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 311d71ae5a4SJacob Faibussowitsch { 3123a40ed3dSBarry Smith PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONLS(snes)); 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3165e76c431SBarry Smith } 31794298653SBarry Smith 318329e7e40SMatthew Knepley /* 31904d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 32004d965bbSLois Curfman McInnes 32104d965bbSLois Curfman McInnes Input Parameters: 32204d965bbSLois Curfman McInnes . SNES - the SNES context 32304d965bbSLois Curfman McInnes . viewer - visualization context 32404d965bbSLois Curfman McInnes 32504d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 32604d965bbSLois Curfman McInnes */ 327d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) 328d71ae5a4SJacob Faibussowitsch { 329ace3abfcSBarry Smith PetscBool iascii; 330a935fc98SLois Curfman McInnes 3313a40ed3dSBarry Smith PetscFunctionBegin; 3329566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3339371c9d4SSatish Balay if (iascii) { } 3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33550f6de3fSJed Brown } 33650f6de3fSJed Brown 33704d965bbSLois Curfman McInnes /* 33804d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 33904d965bbSLois Curfman McInnes 34004d965bbSLois Curfman McInnes Input Parameter: 34104d965bbSLois Curfman McInnes . snes - the SNES context 34204d965bbSLois Curfman McInnes 34304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 34404d965bbSLois Curfman McInnes */ 345d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) 346d71ae5a4SJacob Faibussowitsch { 3473a40ed3dSBarry Smith PetscFunctionBegin; 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3495e42470aSBarry Smith } 3504827ddcaSBarry Smith 351ebe3b25bSBarry Smith /*MC 35204d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 35304d965bbSLois Curfman McInnes 3543c7db156SBarry Smith Options Database Keys: 3555a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3565a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 357f6dfbefdSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`) 35882d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 35982d26c24SPeter 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) 3605a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 36182d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 36282d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 363acbee50cSBarry Smith 364420bcc1bSBarry Smith Level: beginner 365420bcc1bSBarry Smith 366f6dfbefdSBarry Smith Note: 367f6dfbefdSBarry Smith This is the default nonlinear solver in `SNES` 36804d965bbSLois Curfman McInnes 369420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()` 370f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()` 371ebe3b25bSBarry Smith M*/ 372d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 373d71ae5a4SJacob Faibussowitsch { 37404d7464bSBarry Smith SNES_NEWTONLS *neP; 375d8d34be6SBarry Smith SNESLineSearch linesearch; 3765e42470aSBarry Smith 3773a40ed3dSBarry Smith PetscFunctionBegin; 37804d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 37904d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 38004d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 38104d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 38204d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 38304d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 3845e42470aSBarry Smith 385efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 38642f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 387efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 3884fc747eaSLawrence Mitchell 3899566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 39048a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 391d8d34be6SBarry Smith 3924fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3934fc747eaSLawrence Mitchell 3944dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 3955e42470aSBarry Smith snes->data = (void *)neP; 3963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3975e42470aSBarry Smith } 398