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 */ 10f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,PetscReal fnorm,PetscBool *ismin) 118a5d9ee4SBarry Smith { 128a5d9ee4SBarry Smith PetscReal a1; 13ace3abfcSBarry Smith PetscBool hastranspose; 14f6c5ca78SBarry Smith Vec W; 158a5d9ee4SBarry Smith 168a5d9ee4SBarry Smith PetscFunctionBegin; 178a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 185f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(F,&W)); 2036669109SBarry Smith if (hastranspose) { 218a5d9ee4SBarry Smith /* Compute || J^T F|| */ 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,F,W)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(W,NORM_2,&a1)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(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 315f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(W,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(W,NORM_2,&wnorm)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(W,&work)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,W,work)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(F,work,&result)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&work)); 3736669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 385f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&W)); 428a5d9ee4SBarry Smith PetscFunctionReturn(0); 438a5d9ee4SBarry Smith } 448a5d9ee4SBarry Smith 4574637425SBarry Smith /* 465ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4774637425SBarry Smith */ 48f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X) 4974637425SBarry Smith { 5074637425SBarry Smith PetscReal a1,a2; 51ace3abfcSBarry Smith PetscBool hastranspose; 5274637425SBarry Smith 5374637425SBarry Smith PetscFunctionBegin; 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose)); 5574637425SBarry Smith if (hastranspose) { 56f6c5ca78SBarry Smith Vec W1,W2; 57f6c5ca78SBarry Smith 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(F,&W1)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(F,&W2)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,X,W1)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(W1,-1.0,F)); 6274637425SBarry Smith 6374637425SBarry Smith /* Compute || J^T W|| */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,W1,W2)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(W1,NORM_2,&a1)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(W2,NORM_2,&a2)); 6775567043SBarry Smith if (a1 != 0.0) { 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1))); 6985471664SBarry Smith } 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&W1)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&W2)); 7274637425SBarry Smith } 7374637425SBarry Smith PetscFunctionReturn(0); 7474637425SBarry Smith } 7574637425SBarry Smith 7604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7704d965bbSLois Curfman McInnes 7804d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7994b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 8004d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8104d965bbSLois Curfman McInnes respectively. 8204d965bbSLois Curfman McInnes 83fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8404d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8504d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 86fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8704d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8804d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 8904d7464bSBarry Smith we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 904b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 9104d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9204d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9304d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9404d965bbSLois Curfman McInnes for all nonlinear solvers. 9504d965bbSLois Curfman McInnes 9604d965bbSLois Curfman McInnes Another key routine is: 9704d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9804d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9904d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 10004d965bbSLois Curfman McInnes SNESSolve() if necessary. 10104d965bbSLois Curfman McInnes 10204d965bbSLois Curfman McInnes Additional basic routines are: 10304d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10404d965bbSLois Curfman McInnes have actually been used. 10504d965bbSLois Curfman McInnes These are called by application codes via the interface routines 106186905e3SBarry Smith SNESView(). 10704d965bbSLois Curfman McInnes 10804d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10904d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 11004d965bbSLois Curfman McInnes above description applies to these categories also. 11104d965bbSLois Curfman McInnes 11204d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1135e42470aSBarry Smith /* 11404d7464bSBarry Smith SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 11504d965bbSLois Curfman McInnes method with a line search. 1165e76c431SBarry Smith 11704d965bbSLois Curfman McInnes Input Parameters: 11804d965bbSLois Curfman McInnes . snes - the SNES context 1195e76c431SBarry Smith 12004d965bbSLois Curfman McInnes Output Parameter: 121c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12204d965bbSLois Curfman McInnes 12304d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1245e76c431SBarry Smith 1255e76c431SBarry Smith Notes: 1265e76c431SBarry Smith This implements essentially a truncated Newton method with a 1275e76c431SBarry Smith line search. By default a cubic backtracking line search 1285e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1295e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 130393d2d9aSLois Curfman McInnes and Schnabel. 1315e76c431SBarry Smith */ 13204d7464bSBarry Smith PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 1335e76c431SBarry Smith { 134a7cc72afSBarry Smith PetscInt maxits,i,lits; 135422a814eSBarry Smith SNESLineSearchReason lssucceed; 1361936c542SPeter Brune PetscReal fnorm,gnorm,xnorm,ynorm; 137f6c5ca78SBarry Smith Vec Y,X,F; 138f1c6b773SPeter Brune SNESLineSearch linesearch; 139c40d0f55SPeter Brune SNESConvergedReason reason; 1405e76c431SBarry Smith 1413a40ed3dSBarry Smith PetscFunctionBegin; 1422c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 143c579b300SPatrick Farrell 1443d4c4710SBarry Smith snes->numFailures = 0; 1453d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 146184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 147184914b5SBarry Smith 1485e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1495e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 15039e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 151d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1525e76c431SBarry Smith 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1544c49b128SBarry Smith snes->iter = 0; 15575567043SBarry Smith snes->norm = 0.0; 1565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetLineSearch(snes, &linesearch)); 158ce5a860aSPeter Brune 159ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 160efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1615f80ce2aSJacob Faibussowitsch CHKERRQ(SNESApplyNPC(snes,X,NULL,F)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetConvergedReason(snes->npc,&reason)); 163b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 164b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 165b7281c8aSPeter Brune PetscFunctionReturn(0); 166b7281c8aSPeter Brune } 167b7281c8aSPeter Brune 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormBegin(F,NORM_2,&fnorm)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecNormEnd(F,NORM_2,&fnorm)); 170ce5a860aSPeter Brune } else { 171e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1725f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeFunction(snes,X,F)); 1731aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 174c1c75074SPeter Brune } 1751aa26658SKarl Rupp 1765f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(F,NORM_2,&fnorm)); /* fnorm <- ||F|| */ 177422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSAWsTakeAccess((PetscObject)snes)); 1795e42470aSBarry Smith snes->norm = fnorm; 1805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLogConvergenceHistory(snes,fnorm,0)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(SNESMonitor(snes,0,fnorm)); 1833f149594SLisandro Dalcin 18485385478SLisandro Dalcin /* test convergence */ 1855f80ce2aSJacob Faibussowitsch CHKERRQ((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 18606ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 187d034289bSBarry Smith 1885e76c431SBarry Smith for (i=0; i<maxits; i++) { 1895e76c431SBarry Smith 190329e7e40SMatthew Knepley /* Call general purpose update function */ 191e7788613SBarry Smith if (snes->ops->update) { 1925f80ce2aSJacob Faibussowitsch CHKERRQ((*snes->ops->update)(snes, snes->iter)); 193de8cb200SMatthew Knepley } 194329e7e40SMatthew Knepley 195ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 196efd4aadfSBarry Smith if (snes->npc) { 197efd4aadfSBarry Smith if (snes->npcside== PC_RIGHT) { 1985f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetInitialFunction(snes->npc, F)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,snes->vec_rhs,0)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes->npc, snes->vec_rhs, X)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,snes->vec_rhs,0)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetConvergedReason(snes->npc,&reason)); 203c40d0f55SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 204c40d0f55SPeter Brune snes->reason = SNES_DIVERGED_INNER; 205c40d0f55SPeter Brune PetscFunctionReturn(0); 206c40d0f55SPeter Brune } 2075f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetNPCFunction(snes,F,&fnorm)); 208efd4aadfSBarry Smith } else if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2095f80ce2aSJacob Faibussowitsch CHKERRQ(SNESApplyNPC(snes,X,F,F)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetConvergedReason(snes->npc,&reason)); 211b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 212b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 213b7281c8aSPeter Brune PetscFunctionReturn(0); 214b7281c8aSPeter Brune } 215ce5a860aSPeter Brune } 216c40d0f55SPeter Brune } 217c40d0f55SPeter Brune 218ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 2195f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre)); 22007b62357SFande Kong SNESCheckJacobianDomainerror(snes); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(snes->ksp,F,Y)); 223422a814eSBarry Smith SNESCheckKSPSolve(snes); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetIterationNumber(snes->ksp,&lits)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits)); 22674637425SBarry Smith 227b0a32e0cSBarry Smith if (PetscLogPrintInfo) { 2285f80ce2aSJacob Faibussowitsch CHKERRQ(SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y)); 22985471664SBarry Smith } 23074637425SBarry Smith 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 */ 2351936c542SPeter Brune gnorm = fnorm; 2365f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLineSearchApply(linesearch, X, F, &fnorm, Y)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLineSearchGetReason(linesearch, &lssucceed)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)gnorm,(double)fnorm,(double)ynorm,(int)lssucceed)); 2407e49c775SBarry Smith if (snes->reason) break; 241422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 242422a814eSBarry Smith if (lssucceed) { 243c21ba15cSPeter Brune if (snes->stol*xnorm > ynorm) { 244c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 245c21ba15cSPeter Brune PetscFunctionReturn(0); 246c21ba15cSPeter Brune } 24750ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 248ace3abfcSBarry Smith PetscBool ismin; 249647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 2505f80ce2aSJacob Faibussowitsch CHKERRQ(SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,fnorm,&ismin)); 2518a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2523add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2533add74b1SBarry Smith PetscViewer monitor; 2545f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLineSearchGetDefaultMonitor(linesearch,&monitor)); 255*28b400f6SJacob 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]); 25698921bdaSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged due %s.",SNESConvergedReasons[snes->reason]); 2573add74b1SBarry Smith } 2583505fcc1SBarry Smith break; 2593505fcc1SBarry Smith } 26050ffb88aSMatthew Knepley } 26185385478SLisandro Dalcin /* Monitor convergence */ 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSAWsTakeAccess((PetscObject)snes)); 26385385478SLisandro Dalcin snes->iter = i+1; 26485385478SLisandro Dalcin snes->norm = fnorm; 265c1e67a49SFande Kong snes->ynorm = ynorm; 266c1e67a49SFande Kong snes->xnorm = xnorm; 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLogConvergenceHistory(snes,snes->norm,lits)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(SNESMonitor(snes,snes->iter,snes->norm)); 2702a997342SPeter Brune /* Test for convergence */ 2715f80ce2aSJacob Faibussowitsch CHKERRQ((*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP)); 2723f149594SLisandro Dalcin if (snes->reason) break; 27316c95cacSBarry Smith } 27452392280SLois Curfman McInnes if (i == maxits) { 2755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"Maximum number of iterations has been reached: %D\n",maxits)); 27685385478SLisandro Dalcin if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 27752392280SLois Curfman McInnes } 2783a40ed3dSBarry Smith PetscFunctionReturn(0); 2795e76c431SBarry Smith } 28004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 28104d965bbSLois Curfman McInnes /* 28204d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 28304d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 28404d965bbSLois Curfman McInnes 28504d965bbSLois Curfman McInnes Input Parameter: 28604d965bbSLois Curfman McInnes . snes - the SNES context 28704d965bbSLois Curfman McInnes . x - the solution vector 28804d965bbSLois Curfman McInnes 28904d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 29004d965bbSLois Curfman McInnes 29104d965bbSLois Curfman McInnes Notes: 2924b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 29304d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 29404d965bbSLois Curfman McInnes the call to SNESSolve(). 29504d965bbSLois Curfman McInnes */ 29604d7464bSBarry Smith PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 2975e76c431SBarry Smith { 2983a40ed3dSBarry Smith PetscFunctionBegin; 2995f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetUpMatrices(snes)); 300efd4aadfSBarry Smith if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 3013a40ed3dSBarry Smith PetscFunctionReturn(0); 3025e76c431SBarry Smith } 30304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3046b8b9a38SLisandro Dalcin 30504d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONLS(SNES snes) 3066b8b9a38SLisandro Dalcin { 3076b8b9a38SLisandro Dalcin PetscFunctionBegin; 3086b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 3096b8b9a38SLisandro Dalcin } 3106b8b9a38SLisandro Dalcin 31104d965bbSLois Curfman McInnes /* 31204d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 31304d7464bSBarry Smith with SNESCreate_NEWTONLS(). 31404d965bbSLois Curfman McInnes 31504d965bbSLois Curfman McInnes Input Parameter: 31604d965bbSLois Curfman McInnes . snes - the SNES context 31704d965bbSLois Curfman McInnes 318de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 31904d965bbSLois Curfman McInnes */ 32004d7464bSBarry Smith PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 3215e76c431SBarry Smith { 3223a40ed3dSBarry Smith PetscFunctionBegin; 3235f80ce2aSJacob Faibussowitsch CHKERRQ(SNESReset_NEWTONLS(snes)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(snes->data)); 3253a40ed3dSBarry Smith PetscFunctionReturn(0); 3265e76c431SBarry Smith } 32704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 32894298653SBarry Smith 329329e7e40SMatthew Knepley /* 33004d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 33104d965bbSLois Curfman McInnes 33204d965bbSLois Curfman McInnes Input Parameters: 33304d965bbSLois Curfman McInnes . SNES - the SNES context 33404d965bbSLois Curfman McInnes . viewer - visualization context 33504d965bbSLois Curfman McInnes 33604d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 33704d965bbSLois Curfman McInnes */ 33804d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer) 339a935fc98SLois Curfman McInnes { 340ace3abfcSBarry Smith PetscBool iascii; 341a935fc98SLois Curfman McInnes 3423a40ed3dSBarry Smith PetscFunctionBegin; 3435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 34432077d6dSBarry Smith if (iascii) { 34550f6de3fSJed Brown } 34650f6de3fSJed Brown PetscFunctionReturn(0); 34750f6de3fSJed Brown } 34850f6de3fSJed Brown 34904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 35004d965bbSLois Curfman McInnes /* 35104d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 35204d965bbSLois Curfman McInnes 35304d965bbSLois Curfman McInnes Input Parameter: 35404d965bbSLois Curfman McInnes . snes - the SNES context 35504d965bbSLois Curfman McInnes 35604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 35704d965bbSLois Curfman McInnes */ 3584416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(PetscOptionItems *PetscOptionsObject,SNES snes) 3595e42470aSBarry Smith { 3603a40ed3dSBarry Smith PetscFunctionBegin; 3613a40ed3dSBarry Smith PetscFunctionReturn(0); 3625e42470aSBarry Smith } 3634827ddcaSBarry Smith 36404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 365ebe3b25bSBarry Smith /*MC 36604d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 36704d965bbSLois Curfman McInnes 368ebe3b25bSBarry Smith Options Database: 3695a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3705a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 3711fe24845SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (SNESLineSearchSetComputeNorms()) 37282d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 37382d26c24SPeter 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) 3745a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 37582d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 37682d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 377acbee50cSBarry Smith 37895452b02SPatrick Sanan Notes: 37995452b02SPatrick Sanan This is the default nonlinear solver in SNES 38004d965bbSLois Curfman McInnes 381ee3001cbSBarry Smith Level: beginner 382ee3001cbSBarry Smith 38304d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder() 3845a9b6599SPeter Brune SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms() 385ebe3b25bSBarry Smith 386ebe3b25bSBarry Smith M*/ 3878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 3885e42470aSBarry Smith { 38904d7464bSBarry Smith SNES_NEWTONLS *neP; 390d8d34be6SBarry Smith SNESLineSearch linesearch; 3915e42470aSBarry Smith 3923a40ed3dSBarry Smith PetscFunctionBegin; 39304d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 39404d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 39504d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 39604d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 39704d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 39804d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 3995e42470aSBarry Smith 400efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 40142f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 402efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 4034fc747eaSLawrence Mitchell 4045f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetLineSearch(snes, &linesearch)); 405ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 4065f80ce2aSJacob Faibussowitsch CHKERRQ(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT)); 407ec786807SJed Brown } 408d8d34be6SBarry Smith 4094fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 4104fc747eaSLawrence Mitchell 4115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(snes,&neP)); 4125e42470aSBarry Smith snes->data = (void*)neP; 4133a40ed3dSBarry Smith PetscFunctionReturn(0); 4145e42470aSBarry Smith } 415