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)); 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 } 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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; 127*4279555eSSatish Balay PetscReal fnorm, xnorm, ynorm; 128f6c5ca78SBarry Smith Vec Y, X, F; 129f1c6b773SPeter Brune SNESLineSearch linesearch; 130c40d0f55SPeter Brune SNESConvergedReason reason; 131*4279555eSSatish Balay #if defined(PETSC_USE_INFO) 132*4279555eSSatish Balay PetscReal gnorm; 133*4279555eSSatish Balay #endif 1345e76c431SBarry Smith 1353a40ed3dSBarry Smith PetscFunctionBegin; 1360b121fc5SBarry 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); 137c579b300SPatrick Farrell 1383d4c4710SBarry Smith snes->numFailures = 0; 1393d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 140184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 141184914b5SBarry Smith 1425e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1435e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 14439e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 145d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1465e76c431SBarry Smith 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1484c49b128SBarry Smith snes->iter = 0; 14975567043SBarry Smith snes->norm = 0.0; 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1519566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 152ce5a860aSPeter Brune 153ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 154efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1559566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1569566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 157b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 158b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160b7281c8aSPeter Brune } 161b7281c8aSPeter Brune 1629566063dSJacob Faibussowitsch PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 1639566063dSJacob Faibussowitsch PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 164ce5a860aSPeter Brune } else { 165e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1669566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 1671aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 168c1c75074SPeter Brune } 1691aa26658SKarl Rupp 1709566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 171422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 1729566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1735e42470aSBarry Smith snes->norm = fnorm; 1749566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1759566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 1769566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 1773f149594SLisandro Dalcin 17885385478SLisandro Dalcin /* test convergence */ 179dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 1803ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 181d034289bSBarry Smith 1825e76c431SBarry Smith for (i = 0; i < maxits; i++) { 183329e7e40SMatthew Knepley /* Call general purpose update function */ 184dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 185329e7e40SMatthew Knepley 186ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 187efd4aadfSBarry Smith if (snes->npc) { 188efd4aadfSBarry Smith if (snes->npcside == PC_RIGHT) { 1899566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F)); 1909566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1919566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 1929566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 1939566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 194c40d0f55SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 195c40d0f55SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197c40d0f55SPeter Brune } 1989566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 199efd4aadfSBarry Smith } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2009566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, F)); 2019566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 202b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 203b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 205b7281c8aSPeter Brune } 206ce5a860aSPeter Brune } 207c40d0f55SPeter Brune } 208c40d0f55SPeter Brune 209ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 2109566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 21107b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2129566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 2139566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, Y)); 214422a814eSBarry Smith SNESCheckKSPSolve(snes); 2159566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 21663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 21774637425SBarry Smith 2181baa6e33SBarry Smith if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y)); 21974637425SBarry Smith 220*4279555eSSatish Balay #if defined(PETSC_USE_INFO) 221*4279555eSSatish Balay gnorm = fnorm; 222*4279555eSSatish Balay #endif 223ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 224d12b9ff3SPeter Brune X <- X - lambda*Y 225d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 226ea4d3ed3SLois Curfman McInnes */ 2279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 2289566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed)); 2299566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 2309566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed)); 2317e49c775SBarry Smith if (snes->reason) break; 232422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 233422a814eSBarry Smith if (lssucceed) { 234c21ba15cSPeter Brune if (snes->stol * xnorm > ynorm) { 235c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 2363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 237c21ba15cSPeter Brune } 23850ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 239ace3abfcSBarry Smith PetscBool ismin; 240647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2419566063dSJacob Faibussowitsch PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin)); 2428a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2433add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2443add74b1SBarry Smith PetscViewer monitor; 2459566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor)); 24628b400f6SJacob 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]); 247f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]); 2483add74b1SBarry Smith } 2493505fcc1SBarry Smith break; 2503505fcc1SBarry Smith } 25150ffb88aSMatthew Knepley } 25285385478SLisandro Dalcin /* Monitor convergence */ 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 25485385478SLisandro Dalcin snes->iter = i + 1; 25585385478SLisandro Dalcin snes->norm = fnorm; 256c1e67a49SFande Kong snes->ynorm = ynorm; 257c1e67a49SFande Kong snes->xnorm = xnorm; 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2599566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 2609566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2612a997342SPeter Brune /* Test for convergence */ 262dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 2633f149594SLisandro Dalcin if (snes->reason) break; 26416c95cacSBarry Smith } 26552392280SLois Curfman McInnes if (i == maxits) { 26663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 26785385478SLisandro Dalcin if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 26852392280SLois Curfman McInnes } 2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2705e76c431SBarry Smith } 271f6dfbefdSBarry Smith 27204d965bbSLois Curfman McInnes /* 27304d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 27404d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 27504d965bbSLois Curfman McInnes 27604d965bbSLois Curfman McInnes Input Parameter: 27704d965bbSLois Curfman McInnes . snes - the SNES context 27804d965bbSLois Curfman McInnes . x - the solution vector 27904d965bbSLois Curfman McInnes 28004d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 28104d965bbSLois Curfman McInnes 28204d965bbSLois Curfman McInnes */ 283d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 284d71ae5a4SJacob Faibussowitsch { 2853a40ed3dSBarry Smith PetscFunctionBegin; 2869566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 287efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2895e76c431SBarry Smith } 2906b8b9a38SLisandro Dalcin 291d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONLS(SNES snes) 292d71ae5a4SJacob Faibussowitsch { 2936b8b9a38SLisandro Dalcin PetscFunctionBegin; 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2956b8b9a38SLisandro Dalcin } 2966b8b9a38SLisandro Dalcin 29704d965bbSLois Curfman McInnes /* 29804d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 29904d7464bSBarry Smith with SNESCreate_NEWTONLS(). 30004d965bbSLois Curfman McInnes 30104d965bbSLois Curfman McInnes Input Parameter: 30204d965bbSLois Curfman McInnes . snes - the SNES context 30304d965bbSLois Curfman McInnes 304de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 30504d965bbSLois Curfman McInnes */ 306d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 307d71ae5a4SJacob Faibussowitsch { 3083a40ed3dSBarry Smith PetscFunctionBegin; 3099566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONLS(snes)); 3109566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3125e76c431SBarry Smith } 31394298653SBarry Smith 314329e7e40SMatthew Knepley /* 31504d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 31604d965bbSLois Curfman McInnes 31704d965bbSLois Curfman McInnes Input Parameters: 31804d965bbSLois Curfman McInnes . SNES - the SNES context 31904d965bbSLois Curfman McInnes . viewer - visualization context 32004d965bbSLois Curfman McInnes 32104d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 32204d965bbSLois Curfman McInnes */ 323d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) 324d71ae5a4SJacob Faibussowitsch { 325ace3abfcSBarry Smith PetscBool iascii; 326a935fc98SLois Curfman McInnes 3273a40ed3dSBarry Smith PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3299371c9d4SSatish Balay if (iascii) { } 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33150f6de3fSJed Brown } 33250f6de3fSJed Brown 33304d965bbSLois Curfman McInnes /* 33404d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 33504d965bbSLois Curfman McInnes 33604d965bbSLois Curfman McInnes Input Parameter: 33704d965bbSLois Curfman McInnes . snes - the SNES context 33804d965bbSLois Curfman McInnes 33904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 34004d965bbSLois Curfman McInnes */ 341d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) 342d71ae5a4SJacob Faibussowitsch { 3433a40ed3dSBarry Smith PetscFunctionBegin; 3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3455e42470aSBarry Smith } 3464827ddcaSBarry Smith 347ebe3b25bSBarry Smith /*MC 34804d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 34904d965bbSLois Curfman McInnes 3503c7db156SBarry Smith Options Database Keys: 3515a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3525a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 353f6dfbefdSBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`) 35482d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 35582d26c24SPeter 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) 3565a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 35782d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 35882d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 359acbee50cSBarry Smith 360f6dfbefdSBarry Smith Note: 361f6dfbefdSBarry Smith This is the default nonlinear solver in `SNES` 36204d965bbSLois Curfman McInnes 363ee3001cbSBarry Smith Level: beginner 364ee3001cbSBarry Smith 365db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()` 366f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()` 367ebe3b25bSBarry Smith M*/ 368d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 369d71ae5a4SJacob Faibussowitsch { 37004d7464bSBarry Smith SNES_NEWTONLS *neP; 371d8d34be6SBarry Smith SNESLineSearch linesearch; 3725e42470aSBarry Smith 3733a40ed3dSBarry Smith PetscFunctionBegin; 37404d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 37504d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 37604d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 37704d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 37804d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 37904d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 3805e42470aSBarry Smith 381efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 38242f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 383efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 3844fc747eaSLawrence Mitchell 3859566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 38648a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 387d8d34be6SBarry Smith 3884fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 3894fc747eaSLawrence Mitchell 3904dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 3915e42470aSBarry Smith snes->data = (void *)neP; 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3935e42470aSBarry Smith } 394