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 */ 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin) 11d71ae5a4SJacob Faibussowitsch { 128a5d9ee4SBarry Smith PetscReal a1; 13ace3abfcSBarry Smith PetscBool hastranspose; 14f6c5ca78SBarry Smith Vec W; 158a5d9ee4SBarry Smith 168a5d9ee4SBarry Smith PetscFunctionBegin; 178a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 189566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W)); 2036669109SBarry Smith if (hastranspose) { 218a5d9ee4SBarry Smith /* Compute || J^T F|| */ 229566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, F, W)); 239566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &a1)); 249566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm))); 258a5d9ee4SBarry Smith if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE; 2636669109SBarry Smith } else { 2736669109SBarry Smith Vec work; 28ea709b57SSatish Balay PetscScalar result; 2936669109SBarry Smith PetscReal wnorm; 3036669109SBarry Smith 319566063dSJacob Faibussowitsch PetscCall(VecSetRandom(W, NULL)); 329566063dSJacob Faibussowitsch PetscCall(VecNorm(W, NORM_2, &wnorm)); 339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(W, &work)); 349566063dSJacob Faibussowitsch PetscCall(MatMult(A, W, work)); 359566063dSJacob Faibussowitsch PetscCall(VecDot(F, work, &result)); 369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work)); 3736669109SBarry Smith a1 = PetscAbsScalar(result) / (fnorm * wnorm); 389566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1)); 3936669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4036669109SBarry Smith } 419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W)); 428a5d9ee4SBarry Smith PetscFunctionReturn(0); 438a5d9ee4SBarry Smith } 448a5d9ee4SBarry Smith 4574637425SBarry Smith /* 465ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4774637425SBarry Smith */ 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X) 49d71ae5a4SJacob Faibussowitsch { 5074637425SBarry Smith PetscReal a1, a2; 51ace3abfcSBarry Smith PetscBool hastranspose; 5274637425SBarry Smith 5374637425SBarry Smith PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose)); 5574637425SBarry Smith if (hastranspose) { 56f6c5ca78SBarry Smith Vec W1, W2; 57f6c5ca78SBarry Smith 589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W1)); 599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(F, &W2)); 609566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, W1)); 619566063dSJacob Faibussowitsch PetscCall(VecAXPY(W1, -1.0, F)); 6274637425SBarry Smith 6374637425SBarry Smith /* Compute || J^T W|| */ 649566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 659566063dSJacob Faibussowitsch PetscCall(VecNorm(W1, NORM_2, &a1)); 669566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &a2)); 6748a46eb9SPierre 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))); 689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W1)); 699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2)); 7074637425SBarry Smith } 7174637425SBarry Smith PetscFunctionReturn(0); 7274637425SBarry Smith } 7374637425SBarry Smith 74f6dfbefdSBarry Smith /* 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 8604d7464bSBarry Smith we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) 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 109f6dfbefdSBarry Smith */ 1105e42470aSBarry Smith /* 11104d7464bSBarry Smith SNESSolve_NEWTONLS - 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 */ 123d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 124d71ae5a4SJacob Faibussowitsch { 125a7cc72afSBarry Smith PetscInt maxits, i, lits; 126422a814eSBarry Smith SNESLineSearchReason lssucceed; 1271936c542SPeter Brune PetscReal fnorm, gnorm, xnorm, ynorm; 128f6c5ca78SBarry Smith Vec Y, X, F; 129f1c6b773SPeter Brune SNESLineSearch linesearch; 130c40d0f55SPeter Brune SNESConvergedReason reason; 1315e76c431SBarry Smith 1323a40ed3dSBarry Smith PetscFunctionBegin; 1330b121fc5SBarry 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); 134c579b300SPatrick Farrell 1353d4c4710SBarry Smith snes->numFailures = 0; 1363d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 137184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 138184914b5SBarry Smith 1395e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1405e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14139e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 142d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1435e76c431SBarry Smith 1449566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1454c49b128SBarry Smith snes->iter = 0; 14675567043SBarry Smith snes->norm = 0.0; 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1489566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 149ce5a860aSPeter Brune 150ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 151efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1529566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1539566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 154b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 155b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 156b7281c8aSPeter Brune PetscFunctionReturn(0); 157b7281c8aSPeter Brune } 158b7281c8aSPeter Brune 1599566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 1609566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 161ce5a860aSPeter Brune } else { 162e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1639566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 1641aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 165c1c75074SPeter Brune } 1661aa26658SKarl Rupp 1679566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 168422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 1699566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1705e42470aSBarry Smith snes->norm = fnorm; 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1729566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 1739566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 1743f149594SLisandro Dalcin 17585385478SLisandro Dalcin /* test convergence */ 176dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 17706ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 178d034289bSBarry Smith 1795e76c431SBarry Smith for (i = 0; i < maxits; i++) { 180329e7e40SMatthew Knepley /* Call general purpose update function */ 181dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 182329e7e40SMatthew Knepley 183ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 184efd4aadfSBarry Smith if (snes->npc) { 185efd4aadfSBarry Smith if (snes->npcside == PC_RIGHT) { 1869566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F)); 1879566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1889566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 1899566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1909566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 191c40d0f55SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 192c40d0f55SPeter Brune snes->reason = SNES_DIVERGED_INNER; 193c40d0f55SPeter Brune PetscFunctionReturn(0); 194c40d0f55SPeter Brune } 1959566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 196efd4aadfSBarry Smith } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 1979566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, F)); 1989566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 199b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 200b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 201b7281c8aSPeter Brune PetscFunctionReturn(0); 202b7281c8aSPeter Brune } 203ce5a860aSPeter Brune } 204c40d0f55SPeter Brune } 205c40d0f55SPeter Brune 206ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 2079566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 20807b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2099566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2109566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 211422a814eSBarry Smith SNESCheckKSPSolve(snes); 2129566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 21363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 21474637425SBarry Smith 2151baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y)); 21674637425SBarry Smith 217ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 218d12b9ff3SPeter Brune X <- X - lambda*Y 219d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 220ea4d3ed3SLois Curfman McInnes */ 2211936c542SPeter Brune gnorm = fnorm; 2229566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 2239566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed)); 2249566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 2259566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed)); 2267e49c775SBarry Smith if (snes->reason) break; 227422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 228422a814eSBarry Smith if (lssucceed) { 229c21ba15cSPeter Brune if (snes->stol * xnorm > ynorm) { 230c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 231c21ba15cSPeter Brune PetscFunctionReturn(0); 232c21ba15cSPeter Brune } 23350ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 234ace3abfcSBarry Smith PetscBool ismin; 235647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2369566063dSJacob Faibussowitsch PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin)); 2378a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2383add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2393add74b1SBarry Smith PetscViewer monitor; 2409566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 24128b400f6SJacob 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]); 242f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]); 2433add74b1SBarry Smith } 2443505fcc1SBarry Smith break; 2453505fcc1SBarry Smith } 24650ffb88aSMatthew Knepley } 24785385478SLisandro Dalcin /* Monitor convergence */ 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 24985385478SLisandro Dalcin snes->iter = i + 1; 25085385478SLisandro Dalcin snes->norm = fnorm; 251c1e67a49SFande Kong snes->ynorm = ynorm; 252c1e67a49SFande Kong snes->xnorm = xnorm; 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2549566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 2559566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2562a997342SPeter Brune /* Test for convergence */ 257dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 2583f149594SLisandro Dalcin if (snes->reason) break; 25916c95cacSBarry Smith } 26052392280SLois Curfman McInnes if (i == maxits) { 26163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 26285385478SLisandro Dalcin if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 26352392280SLois Curfman McInnes } 2643a40ed3dSBarry Smith PetscFunctionReturn(0); 2655e76c431SBarry Smith } 266f6dfbefdSBarry Smith 26704d965bbSLois Curfman McInnes /* 26804d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 26904d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 27004d965bbSLois Curfman McInnes 27104d965bbSLois Curfman McInnes Input Parameter: 27204d965bbSLois Curfman McInnes . snes - the SNES context 27304d965bbSLois Curfman McInnes . x - the solution vector 27404d965bbSLois Curfman McInnes 27504d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 27604d965bbSLois Curfman McInnes 27704d965bbSLois Curfman McInnes */ 278d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 279d71ae5a4SJacob Faibussowitsch { 2803a40ed3dSBarry Smith PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 282efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 2845e76c431SBarry Smith } 2856b8b9a38SLisandro Dalcin 286d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONLS(SNES snes) 287d71ae5a4SJacob Faibussowitsch { 2886b8b9a38SLisandro Dalcin PetscFunctionBegin; 2896b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 2906b8b9a38SLisandro Dalcin } 2916b8b9a38SLisandro Dalcin 29204d965bbSLois Curfman McInnes /* 29304d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 29404d7464bSBarry Smith with SNESCreate_NEWTONLS(). 29504d965bbSLois Curfman McInnes 29604d965bbSLois Curfman McInnes Input Parameter: 29704d965bbSLois Curfman McInnes . snes - the SNES context 29804d965bbSLois Curfman McInnes 299de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30004d965bbSLois Curfman McInnes */ 301d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 302d71ae5a4SJacob Faibussowitsch { 3033a40ed3dSBarry Smith PetscFunctionBegin; 3049566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONLS(snes)); 3059566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3063a40ed3dSBarry Smith PetscFunctionReturn(0); 3075e76c431SBarry Smith } 30894298653SBarry Smith 309329e7e40SMatthew Knepley /* 31004d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 31104d965bbSLois Curfman McInnes 31204d965bbSLois Curfman McInnes Input Parameters: 31304d965bbSLois Curfman McInnes . SNES - the SNES context 31404d965bbSLois Curfman McInnes . viewer - visualization context 31504d965bbSLois Curfman McInnes 31604d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 31704d965bbSLois Curfman McInnes */ 318d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) 319d71ae5a4SJacob Faibussowitsch { 320ace3abfcSBarry Smith PetscBool iascii; 321a935fc98SLois Curfman McInnes 3223a40ed3dSBarry Smith PetscFunctionBegin; 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3249371c9d4SSatish Balay if (iascii) { } 32550f6de3fSJed Brown PetscFunctionReturn(0); 32650f6de3fSJed Brown } 32750f6de3fSJed Brown 32804d965bbSLois Curfman McInnes /* 32904d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 33004d965bbSLois Curfman McInnes 33104d965bbSLois Curfman McInnes Input Parameter: 33204d965bbSLois Curfman McInnes . snes - the SNES context 33304d965bbSLois Curfman McInnes 33404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 33504d965bbSLois Curfman McInnes */ 336d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) 337d71ae5a4SJacob Faibussowitsch { 3383a40ed3dSBarry Smith PetscFunctionBegin; 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 3405e42470aSBarry Smith } 3414827ddcaSBarry Smith 342ebe3b25bSBarry Smith /*MC 34304d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 34404d965bbSLois Curfman McInnes 345*3c7db156SBarry Smith Options Database Keys: 3465a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3475a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 348f6dfbefdSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`) 34982d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 35082d26c24SPeter 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) 3515a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 35282d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 35382d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 354acbee50cSBarry Smith 355f6dfbefdSBarry Smith Note: 356f6dfbefdSBarry Smith This is the default nonlinear solver in `SNES` 35704d965bbSLois Curfman McInnes 358ee3001cbSBarry Smith Level: beginner 359ee3001cbSBarry Smith 360db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()` 361f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()` 362ebe3b25bSBarry Smith M*/ 363d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 364d71ae5a4SJacob Faibussowitsch { 36504d7464bSBarry Smith SNES_NEWTONLS *neP; 366d8d34be6SBarry Smith SNESLineSearch linesearch; 3675e42470aSBarry Smith 3683a40ed3dSBarry Smith PetscFunctionBegin; 36904d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 37004d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 37104d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 37204d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 37304d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 37404d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 3755e42470aSBarry Smith 376efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 37742f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 378efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 3794fc747eaSLawrence Mitchell 3809566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 38148a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 382d8d34be6SBarry Smith 3834fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3844fc747eaSLawrence Mitchell 3854dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 3865e42470aSBarry Smith snes->data = (void *)neP; 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 3885e42470aSBarry Smith } 389