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; 13dfbe8321SBarry Smith PetscErrorCode ierr; 14ace3abfcSBarry Smith PetscBool hastranspose; 15f6c5ca78SBarry Smith Vec W; 168a5d9ee4SBarry Smith 178a5d9ee4SBarry Smith PetscFunctionBegin; 188a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 1936669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 20f6c5ca78SBarry Smith ierr = VecDuplicate(F,&W);CHKERRQ(ierr); 2136669109SBarry Smith if (hastranspose) { 228a5d9ee4SBarry Smith /* Compute || J^T F|| */ 238a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 248a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 257d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 268a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2736669109SBarry Smith } else { 2836669109SBarry Smith Vec work; 29ea709b57SSatish Balay PetscScalar result; 3036669109SBarry Smith PetscReal wnorm; 3136669109SBarry Smith 320298fd71SBarry Smith ierr = VecSetRandom(W,NULL);CHKERRQ(ierr); 3336669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3536669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 376bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 3836669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 397d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 4036669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4136669109SBarry Smith } 42f6c5ca78SBarry Smith ierr = VecDestroy(&W);CHKERRQ(ierr); 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 475ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4874637425SBarry Smith */ 49f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X) 5074637425SBarry Smith { 5174637425SBarry Smith PetscReal a1,a2; 52dfbe8321SBarry Smith PetscErrorCode ierr; 53ace3abfcSBarry Smith PetscBool hastranspose; 5474637425SBarry Smith 5574637425SBarry Smith PetscFunctionBegin; 5674637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 5774637425SBarry Smith if (hastranspose) { 58f6c5ca78SBarry Smith Vec W1,W2; 59f6c5ca78SBarry Smith 60f6c5ca78SBarry Smith ierr = VecDuplicate(F,&W1);CHKERRQ(ierr); 61f6c5ca78SBarry Smith ierr = VecDuplicate(F,&W2);CHKERRQ(ierr); 6274637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6379f36460SBarry Smith ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 6474637425SBarry Smith 6574637425SBarry Smith /* Compute || J^T W|| */ 6674637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6774637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6874637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 6975567043SBarry Smith if (a1 != 0.0) { 707d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 7185471664SBarry Smith } 72f6c5ca78SBarry Smith ierr = VecDestroy(&W1);CHKERRQ(ierr); 73f6c5ca78SBarry Smith ierr = VecDestroy(&W2);CHKERRQ(ierr); 7474637425SBarry Smith } 7574637425SBarry Smith PetscFunctionReturn(0); 7674637425SBarry Smith } 7774637425SBarry Smith 7804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7904d965bbSLois Curfman McInnes 8004d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 8194b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 8204d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8304d965bbSLois Curfman McInnes respectively. 8404d965bbSLois Curfman McInnes 85fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8604d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8704d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 88fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8904d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 9004d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 9104d7464bSBarry Smith we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 924b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 9304d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9404d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9504d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9604d965bbSLois Curfman McInnes for all nonlinear solvers. 9704d965bbSLois Curfman McInnes 9804d965bbSLois Curfman McInnes Another key routine is: 9904d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 10004d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 10104d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 10204d965bbSLois Curfman McInnes SNESSolve() if necessary. 10304d965bbSLois Curfman McInnes 10404d965bbSLois Curfman McInnes Additional basic routines are: 10504d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10604d965bbSLois Curfman McInnes have actually been used. 10704d965bbSLois Curfman McInnes These are called by application codes via the interface routines 108186905e3SBarry Smith SNESView(). 10904d965bbSLois Curfman McInnes 11004d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 11104d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 11204d965bbSLois Curfman McInnes above description applies to these categories also. 11304d965bbSLois Curfman McInnes 11404d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1155e42470aSBarry Smith /* 11604d7464bSBarry Smith SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 11704d965bbSLois Curfman McInnes method with a line search. 1185e76c431SBarry Smith 11904d965bbSLois Curfman McInnes Input Parameters: 12004d965bbSLois Curfman McInnes . snes - the SNES context 1215e76c431SBarry Smith 12204d965bbSLois Curfman McInnes Output Parameter: 123c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12404d965bbSLois Curfman McInnes 12504d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1265e76c431SBarry Smith 1275e76c431SBarry Smith Notes: 1285e76c431SBarry Smith This implements essentially a truncated Newton method with a 1295e76c431SBarry Smith line search. By default a cubic backtracking line search 1305e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1315e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 132393d2d9aSLois Curfman McInnes and Schnabel. 1335e76c431SBarry Smith */ 13404d7464bSBarry Smith PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 1355e76c431SBarry Smith { 1366849ba73SBarry Smith PetscErrorCode ierr; 137a7cc72afSBarry Smith PetscInt maxits,i,lits; 138422a814eSBarry Smith SNESLineSearchReason lssucceed; 1391936c542SPeter Brune PetscReal fnorm,gnorm,xnorm,ynorm; 140f6c5ca78SBarry Smith Vec Y,X,F; 141f1c6b773SPeter Brune SNESLineSearch linesearch; 142c40d0f55SPeter Brune SNESConvergedReason reason; 1435e76c431SBarry Smith 1443a40ed3dSBarry Smith PetscFunctionBegin; 145*2c71b3e2SJacob 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); 146c579b300SPatrick Farrell 1473d4c4710SBarry Smith snes->numFailures = 0; 1483d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 149184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 150184914b5SBarry Smith 1515e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1525e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 15339e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 154d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1555e76c431SBarry Smith 156e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1574c49b128SBarry Smith snes->iter = 0; 15875567043SBarry Smith snes->norm = 0.0; 159e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1607601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 161ce5a860aSPeter Brune 162ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 163efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 164be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 165efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 166b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 167b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 168b7281c8aSPeter Brune PetscFunctionReturn(0); 169b7281c8aSPeter Brune } 170b7281c8aSPeter Brune 171ce5a860aSPeter Brune ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); 172ce5a860aSPeter Brune ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 173ce5a860aSPeter Brune } else { 174e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 17519717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1761aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 177c1c75074SPeter Brune } 1781aa26658SKarl Rupp 179c1c75074SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 180422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 181e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1825e42470aSBarry Smith snes->norm = fnorm; 183e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 184a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 1857a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1863f149594SLisandro Dalcin 18785385478SLisandro Dalcin /* test convergence */ 18885385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 18906ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 190d034289bSBarry Smith 1915e76c431SBarry Smith for (i=0; i<maxits; i++) { 1925e76c431SBarry Smith 193329e7e40SMatthew Knepley /* Call general purpose update function */ 194e7788613SBarry Smith if (snes->ops->update) { 195e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 196de8cb200SMatthew Knepley } 197329e7e40SMatthew Knepley 198ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 199efd4aadfSBarry Smith if (snes->npc) { 200efd4aadfSBarry Smith if (snes->npcside== PC_RIGHT) { 201efd4aadfSBarry Smith ierr = SNESSetInitialFunction(snes->npc, F);CHKERRQ(ierr); 202efd4aadfSBarry Smith ierr = PetscLogEventBegin(SNES_NPCSolve,snes->npc,X,snes->vec_rhs,0);CHKERRQ(ierr); 203efd4aadfSBarry Smith ierr = SNESSolve(snes->npc, snes->vec_rhs, X);CHKERRQ(ierr); 204efd4aadfSBarry Smith ierr = PetscLogEventEnd(SNES_NPCSolve,snes->npc,X,snes->vec_rhs,0);CHKERRQ(ierr); 205efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 206c40d0f55SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 207c40d0f55SPeter Brune snes->reason = SNES_DIVERGED_INNER; 208c40d0f55SPeter Brune PetscFunctionReturn(0); 209c40d0f55SPeter Brune } 210be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 211efd4aadfSBarry Smith } else if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 212be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,F);CHKERRQ(ierr); 213efd4aadfSBarry Smith ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr); 214b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 215b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 216b7281c8aSPeter Brune PetscFunctionReturn(0); 217b7281c8aSPeter Brune } 218ce5a860aSPeter Brune } 219c40d0f55SPeter Brune } 220c40d0f55SPeter Brune 221ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 222d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 22307b62357SFande Kong SNESCheckJacobianDomainerror(snes); 22423ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 225d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 226422a814eSBarry Smith SNESCheckKSPSolve(snes); 2273d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2287d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 22974637425SBarry Smith 230b0a32e0cSBarry Smith if (PetscLogPrintInfo) { 231f6c5ca78SBarry Smith ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y);CHKERRQ(ierr); 23285471664SBarry Smith } 23374637425SBarry Smith 234ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 235d12b9ff3SPeter Brune X <- X - lambda*Y 236d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 237ea4d3ed3SLois Curfman McInnes */ 2381936c542SPeter Brune gnorm = fnorm; 239f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 240422a814eSBarry Smith ierr = SNESLineSearchGetReason(linesearch, &lssucceed);CHKERRQ(ierr); 241f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 2427d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)gnorm,(double)fnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 2437e49c775SBarry Smith if (snes->reason) break; 244422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 245422a814eSBarry Smith if (lssucceed) { 246c21ba15cSPeter Brune if (snes->stol*xnorm > ynorm) { 247c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 248c21ba15cSPeter Brune PetscFunctionReturn(0); 249c21ba15cSPeter Brune } 25050ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 251ace3abfcSBarry Smith PetscBool ismin; 252647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 253f6c5ca78SBarry Smith ierr = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,fnorm,&ismin);CHKERRQ(ierr); 2548a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2553add74b1SBarry Smith if (snes->errorifnotconverged && snes->reason) { 2563add74b1SBarry Smith PetscViewer monitor; 2573add74b1SBarry Smith ierr = SNESLineSearchGetDefaultMonitor(linesearch,&monitor);CHKERRQ(ierr); 258*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!monitor,PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged due to %s. Suggest running with -snes_linesearch_monitor",SNESConvergedReasons[snes->reason]); 25998921bdaSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged due %s.",SNESConvergedReasons[snes->reason]); 2603add74b1SBarry Smith } 2613505fcc1SBarry Smith break; 2623505fcc1SBarry Smith } 26350ffb88aSMatthew Knepley } 26485385478SLisandro Dalcin /* Monitor convergence */ 265e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 26685385478SLisandro Dalcin snes->iter = i+1; 26785385478SLisandro Dalcin snes->norm = fnorm; 268c1e67a49SFande Kong snes->ynorm = ynorm; 269c1e67a49SFande Kong snes->xnorm = xnorm; 270e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 271a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 2727a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 2732a997342SPeter Brune /* Test for convergence */ 274e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2753f149594SLisandro Dalcin if (snes->reason) break; 27616c95cacSBarry Smith } 27752392280SLois Curfman McInnes if (i == maxits) { 2787d3de750SJacob Faibussowitsch ierr = PetscInfo(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 27985385478SLisandro Dalcin if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 28052392280SLois Curfman McInnes } 2813a40ed3dSBarry Smith PetscFunctionReturn(0); 2825e76c431SBarry Smith } 28304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 28404d965bbSLois Curfman McInnes /* 28504d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 28604d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 28704d965bbSLois Curfman McInnes 28804d965bbSLois Curfman McInnes Input Parameter: 28904d965bbSLois Curfman McInnes . snes - the SNES context 29004d965bbSLois Curfman McInnes . x - the solution vector 29104d965bbSLois Curfman McInnes 29204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 29304d965bbSLois Curfman McInnes 29404d965bbSLois Curfman McInnes Notes: 2954b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 29604d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 29704d965bbSLois Curfman McInnes the call to SNESSolve(). 29804d965bbSLois Curfman McInnes */ 29904d7464bSBarry Smith PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 3005e76c431SBarry Smith { 301dfbe8321SBarry Smith PetscErrorCode ierr; 3023a40ed3dSBarry Smith 3033a40ed3dSBarry Smith PetscFunctionBegin; 3046cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 305efd4aadfSBarry Smith if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 3063a40ed3dSBarry Smith PetscFunctionReturn(0); 3075e76c431SBarry Smith } 30804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3096b8b9a38SLisandro Dalcin 31004d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONLS(SNES snes) 3116b8b9a38SLisandro Dalcin { 3126b8b9a38SLisandro Dalcin PetscFunctionBegin; 3136b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 3146b8b9a38SLisandro Dalcin } 3156b8b9a38SLisandro Dalcin 31604d965bbSLois Curfman McInnes /* 31704d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 31804d7464bSBarry Smith with SNESCreate_NEWTONLS(). 31904d965bbSLois Curfman McInnes 32004d965bbSLois Curfman McInnes Input Parameter: 32104d965bbSLois Curfman McInnes . snes - the SNES context 32204d965bbSLois Curfman McInnes 323de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 32404d965bbSLois Curfman McInnes */ 32504d7464bSBarry Smith PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 3265e76c431SBarry Smith { 327dfbe8321SBarry Smith PetscErrorCode ierr; 3283a40ed3dSBarry Smith 3293a40ed3dSBarry Smith PetscFunctionBegin; 33004d7464bSBarry Smith ierr = SNESReset_NEWTONLS(snes);CHKERRQ(ierr); 3315c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 3323a40ed3dSBarry Smith PetscFunctionReturn(0); 3335e76c431SBarry Smith } 33404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 33594298653SBarry Smith 336329e7e40SMatthew Knepley /* 33704d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 33804d965bbSLois Curfman McInnes 33904d965bbSLois Curfman McInnes Input Parameters: 34004d965bbSLois Curfman McInnes . SNES - the SNES context 34104d965bbSLois Curfman McInnes . viewer - visualization context 34204d965bbSLois Curfman McInnes 34304d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 34404d965bbSLois Curfman McInnes */ 34504d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer) 346a935fc98SLois Curfman McInnes { 347dfbe8321SBarry Smith PetscErrorCode ierr; 348ace3abfcSBarry Smith PetscBool iascii; 349a935fc98SLois Curfman McInnes 3503a40ed3dSBarry Smith PetscFunctionBegin; 351251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 35232077d6dSBarry Smith if (iascii) { 35350f6de3fSJed Brown } 35450f6de3fSJed Brown PetscFunctionReturn(0); 35550f6de3fSJed Brown } 35650f6de3fSJed Brown 35704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 35804d965bbSLois Curfman McInnes /* 35904d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 36004d965bbSLois Curfman McInnes 36104d965bbSLois Curfman McInnes Input Parameter: 36204d965bbSLois Curfman McInnes . snes - the SNES context 36304d965bbSLois Curfman McInnes 36404d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 36504d965bbSLois Curfman McInnes */ 3664416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(PetscOptionItems *PetscOptionsObject,SNES snes) 3675e42470aSBarry Smith { 3683a40ed3dSBarry Smith PetscFunctionBegin; 3693a40ed3dSBarry Smith PetscFunctionReturn(0); 3705e42470aSBarry Smith } 3714827ddcaSBarry Smith 37204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 373ebe3b25bSBarry Smith /*MC 37404d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 37504d965bbSLois Curfman McInnes 376ebe3b25bSBarry Smith Options Database: 3775a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3785a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 3791fe24845SBarry Smith . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (SNESLineSearchSetComputeNorms()) 38082d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 38182d26c24SPeter 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) 3825a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 38382d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 38482d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 385acbee50cSBarry Smith 38695452b02SPatrick Sanan Notes: 38795452b02SPatrick Sanan This is the default nonlinear solver in SNES 38804d965bbSLois Curfman McInnes 389ee3001cbSBarry Smith Level: beginner 390ee3001cbSBarry Smith 39104d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder() 3925a9b6599SPeter Brune SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms() 393ebe3b25bSBarry Smith 394ebe3b25bSBarry Smith M*/ 3958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 3965e42470aSBarry Smith { 397dfbe8321SBarry Smith PetscErrorCode ierr; 39804d7464bSBarry Smith SNES_NEWTONLS *neP; 399d8d34be6SBarry Smith SNESLineSearch linesearch; 4005e42470aSBarry Smith 4013a40ed3dSBarry Smith PetscFunctionBegin; 40204d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 40304d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 40404d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 40504d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 40604d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 40704d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 4085e42470aSBarry Smith 409efd4aadfSBarry Smith snes->npcside = PC_RIGHT; 41042f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 411efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 4124fc747eaSLawrence Mitchell 413d8d34be6SBarry Smith ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 414ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 415d8d34be6SBarry Smith ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 416ec786807SJed Brown } 417d8d34be6SBarry Smith 4184fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 4194fc747eaSLawrence Mitchell 420b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 4215e42470aSBarry Smith snes->data = (void*)neP; 4223a40ed3dSBarry Smith PetscFunctionReturn(0); 4235e42470aSBarry Smith } 424