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 */ 104a2ae208SSatish Balay #undef __FUNCT__ 1104d7464bSBarry Smith #define __FUNCT__ "SNESNEWTONLSCheckLocalMin_Private" 12f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,PetscReal fnorm,PetscBool *ismin) 138a5d9ee4SBarry Smith { 148a5d9ee4SBarry Smith PetscReal a1; 15dfbe8321SBarry Smith PetscErrorCode ierr; 16ace3abfcSBarry Smith PetscBool hastranspose; 17f6c5ca78SBarry Smith Vec W; 188a5d9ee4SBarry Smith 198a5d9ee4SBarry Smith PetscFunctionBegin; 208a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2136669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 22f6c5ca78SBarry Smith ierr = VecDuplicate(F,&W);CHKERRQ(ierr); 2336669109SBarry Smith if (hastranspose) { 248a5d9ee4SBarry Smith /* Compute || J^T F|| */ 258a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 268a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 278f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr); 288a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2936669109SBarry Smith } else { 3036669109SBarry Smith Vec work; 31ea709b57SSatish Balay PetscScalar result; 3236669109SBarry Smith PetscReal wnorm; 3336669109SBarry Smith 340298fd71SBarry Smith ierr = VecSetRandom(W,NULL);CHKERRQ(ierr); 3536669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3636669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3736669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3836669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 396bf464f9SBarry Smith ierr = VecDestroy(&work);CHKERRQ(ierr); 4036669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 418f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr); 4236669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4336669109SBarry Smith } 44f6c5ca78SBarry Smith ierr = VecDestroy(&W);CHKERRQ(ierr); 458a5d9ee4SBarry Smith PetscFunctionReturn(0); 468a5d9ee4SBarry Smith } 478a5d9ee4SBarry Smith 4874637425SBarry Smith /* 495ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 5074637425SBarry Smith */ 514a2ae208SSatish Balay #undef __FUNCT__ 5204d7464bSBarry Smith #define __FUNCT__ "SNESNEWTONLSCheckResidual_Private" 53f6c5ca78SBarry Smith static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X) 5474637425SBarry Smith { 5574637425SBarry Smith PetscReal a1,a2; 56dfbe8321SBarry Smith PetscErrorCode ierr; 57ace3abfcSBarry Smith PetscBool hastranspose; 5874637425SBarry Smith 5974637425SBarry Smith PetscFunctionBegin; 6074637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 6174637425SBarry Smith if (hastranspose) { 62f6c5ca78SBarry Smith Vec W1,W2; 63f6c5ca78SBarry Smith 64f6c5ca78SBarry Smith ierr = VecDuplicate(F,&W1);CHKERRQ(ierr); 65f6c5ca78SBarry Smith ierr = VecDuplicate(F,&W2);CHKERRQ(ierr); 6674637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6779f36460SBarry Smith ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 6874637425SBarry Smith 6974637425SBarry Smith /* Compute || J^T W|| */ 7074637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 7174637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 7274637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 7375567043SBarry Smith if (a1 != 0.0) { 748f1a2a5eSBarry Smith ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr); 7585471664SBarry Smith } 76f6c5ca78SBarry Smith ierr = VecDestroy(&W1);CHKERRQ(ierr); 77f6c5ca78SBarry Smith ierr = VecDestroy(&W2);CHKERRQ(ierr); 7874637425SBarry Smith } 7974637425SBarry Smith PetscFunctionReturn(0); 8074637425SBarry Smith } 8174637425SBarry Smith 8204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 8304d965bbSLois Curfman McInnes 8404d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 8594b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 8604d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8704d965bbSLois Curfman McInnes respectively. 8804d965bbSLois Curfman McInnes 89fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 9004d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 9104d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 92fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 9304d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 9404d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 9504d7464bSBarry Smith we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving 964b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 9704d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9804d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9904d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 10004d965bbSLois Curfman McInnes for all nonlinear solvers. 10104d965bbSLois Curfman McInnes 10204d965bbSLois Curfman McInnes Another key routine is: 10304d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 10404d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 10504d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 10604d965bbSLois Curfman McInnes SNESSolve() if necessary. 10704d965bbSLois Curfman McInnes 10804d965bbSLois Curfman McInnes Additional basic routines are: 10904d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 11004d965bbSLois Curfman McInnes have actually been used. 11104d965bbSLois Curfman McInnes These are called by application codes via the interface routines 112186905e3SBarry Smith SNESView(). 11304d965bbSLois Curfman McInnes 11404d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 11504d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 11604d965bbSLois Curfman McInnes above description applies to these categories also. 11704d965bbSLois Curfman McInnes 11804d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1195e42470aSBarry Smith /* 12004d7464bSBarry Smith SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton 12104d965bbSLois Curfman McInnes method with a line search. 1225e76c431SBarry Smith 12304d965bbSLois Curfman McInnes Input Parameters: 12404d965bbSLois Curfman McInnes . snes - the SNES context 1255e76c431SBarry Smith 12604d965bbSLois Curfman McInnes Output Parameter: 127c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12804d965bbSLois Curfman McInnes 12904d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1305e76c431SBarry Smith 1315e76c431SBarry Smith Notes: 1325e76c431SBarry Smith This implements essentially a truncated Newton method with a 1335e76c431SBarry Smith line search. By default a cubic backtracking line search 1345e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1355e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 136393d2d9aSLois Curfman McInnes and Schnabel. 1375e76c431SBarry Smith */ 1384a2ae208SSatish Balay #undef __FUNCT__ 13904d7464bSBarry Smith #define __FUNCT__ "SNESSolve_NEWTONLS" 14004d7464bSBarry Smith PetscErrorCode SNESSolve_NEWTONLS(SNES snes) 1415e76c431SBarry Smith { 1426849ba73SBarry Smith PetscErrorCode ierr; 143a7cc72afSBarry Smith PetscInt maxits,i,lits; 144*422a814eSBarry Smith SNESLineSearchReason lssucceed; 1451936c542SPeter Brune PetscReal fnorm,gnorm,xnorm,ynorm; 146f6c5ca78SBarry Smith Vec Y,X,F; 147f1c6b773SPeter Brune SNESLineSearch linesearch; 148c40d0f55SPeter Brune SNESConvergedReason reason; 1495e76c431SBarry Smith 1503a40ed3dSBarry Smith PetscFunctionBegin; 151c579b300SPatrick Farrell 152c579b300SPatrick Farrell if (snes->xl || snes->xu || snes->ops->computevariablebounds) { 153c579b300SPatrick Farrell SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 154c579b300SPatrick Farrell } 155c579b300SPatrick Farrell 1563d4c4710SBarry Smith snes->numFailures = 0; 1573d4c4710SBarry Smith snes->numLinearSolveFailures = 0; 158184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 159184914b5SBarry Smith 1605e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1615e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 16239e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 163d12b9ff3SPeter Brune Y = snes->vec_sol_update; /* newton step */ 1645e76c431SBarry Smith 165e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1664c49b128SBarry Smith snes->iter = 0; 16775567043SBarry Smith snes->norm = 0.0; 168e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1697601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 170ce5a860aSPeter Brune 171ce5a860aSPeter Brune /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */ 172ce5a860aSPeter Brune if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 173be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr); 174b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 175b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 176b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 177b7281c8aSPeter Brune PetscFunctionReturn(0); 178b7281c8aSPeter Brune } 179b7281c8aSPeter Brune 180ce5a860aSPeter Brune ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); 181ce5a860aSPeter Brune ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 182ce5a860aSPeter Brune } else { 183e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 18419717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1851aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 186c1c75074SPeter Brune } 1871aa26658SKarl Rupp 188c1c75074SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 189*422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 190e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1915e42470aSBarry Smith snes->norm = fnorm; 192e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 193a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 1947a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1953f149594SLisandro Dalcin 19685385478SLisandro Dalcin /* test convergence */ 19785385478SLisandro Dalcin ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 19806ee9f85SBarry Smith if (snes->reason) PetscFunctionReturn(0); 199d034289bSBarry Smith 2005e76c431SBarry Smith for (i=0; i<maxits; i++) { 2015e76c431SBarry Smith 202329e7e40SMatthew Knepley /* Call general purpose update function */ 203e7788613SBarry Smith if (snes->ops->update) { 204e7788613SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 205de8cb200SMatthew Knepley } 206329e7e40SMatthew Knepley 207ce5a860aSPeter Brune /* apply the nonlinear preconditioner */ 208ce5a860aSPeter Brune if (snes->pc) { 209ce5a860aSPeter Brune if (snes->pcside == PC_RIGHT) { 2103024506fSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 21163e7833aSPeter Brune ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr); 212c40d0f55SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, X);CHKERRQ(ierr); 21363e7833aSPeter Brune ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);CHKERRQ(ierr); 214c40d0f55SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 215c40d0f55SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 216c40d0f55SPeter Brune snes->reason = SNES_DIVERGED_INNER; 217c40d0f55SPeter Brune PetscFunctionReturn(0); 218c40d0f55SPeter Brune } 219be95d8f1SBarry Smith ierr = SNESGetNPCFunction(snes,F,&fnorm);CHKERRQ(ierr); 220ce5a860aSPeter Brune } else if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 221be95d8f1SBarry Smith ierr = SNESApplyNPC(snes,X,F,F);CHKERRQ(ierr); 222b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 223b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 224b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 225b7281c8aSPeter Brune PetscFunctionReturn(0); 226b7281c8aSPeter Brune } 227ce5a860aSPeter Brune } 228c40d0f55SPeter Brune } 229c40d0f55SPeter Brune 230ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 231d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 23223ee1639SBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 233d4211eb9SBarry Smith ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 234*422a814eSBarry Smith SNESCheckKSPSolve(snes); 2353d4c4710SBarry Smith ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 2363f149594SLisandro Dalcin snes->linear_its += lits; 2373f149594SLisandro Dalcin ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 23874637425SBarry Smith 239b0a32e0cSBarry Smith if (PetscLogPrintInfo) { 240f6c5ca78SBarry Smith ierr = SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y);CHKERRQ(ierr); 24185471664SBarry Smith } 24274637425SBarry Smith 243ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 244d12b9ff3SPeter Brune X <- X - lambda*Y 245d12b9ff3SPeter Brune and evaluate F = function(X) (depends on the line search). 246ea4d3ed3SLois Curfman McInnes */ 2471936c542SPeter Brune gnorm = fnorm; 248f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 249*422a814eSBarry Smith ierr = SNESLineSearchGetReason(linesearch, &lssucceed);CHKERRQ(ierr); 250f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 2511936c542SPeter Brune ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)gnorm,(double)fnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr); 2524a93e4ddSBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 253*422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 254*422a814eSBarry Smith if (lssucceed) { 255c21ba15cSPeter Brune if (snes->stol*xnorm > ynorm) { 256c21ba15cSPeter Brune snes->reason = SNES_CONVERGED_SNORM_RELATIVE; 257c21ba15cSPeter Brune PetscFunctionReturn(0); 258c21ba15cSPeter Brune } 25950ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 260ace3abfcSBarry Smith PetscBool ismin; 261647a2e1fSBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 262f6c5ca78SBarry Smith ierr = SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,fnorm,&ismin);CHKERRQ(ierr); 2638a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2643505fcc1SBarry Smith break; 2653505fcc1SBarry Smith } 26650ffb88aSMatthew Knepley } 26785385478SLisandro Dalcin /* Monitor convergence */ 268e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 26985385478SLisandro Dalcin snes->iter = i+1; 27085385478SLisandro Dalcin snes->norm = fnorm; 271e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 272a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 2737a03ce2fSLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 2742a997342SPeter Brune /* Test for convergence */ 275e7788613SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2763f149594SLisandro Dalcin if (snes->reason) break; 27716c95cacSBarry Smith } 27852392280SLois Curfman McInnes if (i == maxits) { 279ae15b995SBarry Smith ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 28085385478SLisandro Dalcin if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 28152392280SLois Curfman McInnes } 2823a40ed3dSBarry Smith PetscFunctionReturn(0); 2835e76c431SBarry Smith } 28404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 28504d965bbSLois Curfman McInnes /* 28604d7464bSBarry Smith SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use 28704d7464bSBarry Smith of the SNESNEWTONLS nonlinear solver. 28804d965bbSLois Curfman McInnes 28904d965bbSLois Curfman McInnes Input Parameter: 29004d965bbSLois Curfman McInnes . snes - the SNES context 29104d965bbSLois Curfman McInnes . x - the solution vector 29204d965bbSLois Curfman McInnes 29304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 29404d965bbSLois Curfman McInnes 29504d965bbSLois Curfman McInnes Notes: 2964b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 29704d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 29804d965bbSLois Curfman McInnes the call to SNESSolve(). 29904d965bbSLois Curfman McInnes */ 3004a2ae208SSatish Balay #undef __FUNCT__ 30104d7464bSBarry Smith #define __FUNCT__ "SNESSetUp_NEWTONLS" 30204d7464bSBarry Smith PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) 3035e76c431SBarry Smith { 304dfbe8321SBarry Smith PetscErrorCode ierr; 3053a40ed3dSBarry Smith 3063a40ed3dSBarry Smith PetscFunctionBegin; 3076cab3a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 3086c67d002SPeter Brune if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 3093a40ed3dSBarry Smith PetscFunctionReturn(0); 3105e76c431SBarry Smith } 31104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3126b8b9a38SLisandro Dalcin 3136b8b9a38SLisandro Dalcin #undef __FUNCT__ 31404d7464bSBarry Smith #define __FUNCT__ "SNESReset_NEWTONLS" 31504d7464bSBarry Smith PetscErrorCode SNESReset_NEWTONLS(SNES snes) 3166b8b9a38SLisandro Dalcin { 3176b8b9a38SLisandro Dalcin PetscFunctionBegin; 3186b8b9a38SLisandro Dalcin PetscFunctionReturn(0); 3196b8b9a38SLisandro Dalcin } 3206b8b9a38SLisandro Dalcin 32104d965bbSLois Curfman McInnes /* 32204d7464bSBarry Smith SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created 32304d7464bSBarry Smith with SNESCreate_NEWTONLS(). 32404d965bbSLois Curfman McInnes 32504d965bbSLois Curfman McInnes Input Parameter: 32604d965bbSLois Curfman McInnes . snes - the SNES context 32704d965bbSLois Curfman McInnes 328de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 32904d965bbSLois Curfman McInnes */ 3304a2ae208SSatish Balay #undef __FUNCT__ 33104d7464bSBarry Smith #define __FUNCT__ "SNESDestroy_NEWTONLS" 33204d7464bSBarry Smith PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) 3335e76c431SBarry Smith { 334dfbe8321SBarry Smith PetscErrorCode ierr; 3353a40ed3dSBarry Smith 3363a40ed3dSBarry Smith PetscFunctionBegin; 33704d7464bSBarry Smith ierr = SNESReset_NEWTONLS(snes);CHKERRQ(ierr); 3385c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 3405e76c431SBarry Smith } 34104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 34294298653SBarry Smith 343329e7e40SMatthew Knepley /* 34404d7464bSBarry Smith SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure. 34504d965bbSLois Curfman McInnes 34604d965bbSLois Curfman McInnes Input Parameters: 34704d965bbSLois Curfman McInnes . SNES - the SNES context 34804d965bbSLois Curfman McInnes . viewer - visualization context 34904d965bbSLois Curfman McInnes 35004d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 35104d965bbSLois Curfman McInnes */ 3524a2ae208SSatish Balay #undef __FUNCT__ 35304d7464bSBarry Smith #define __FUNCT__ "SNESView_NEWTONLS" 35404d7464bSBarry Smith static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer) 355a935fc98SLois Curfman McInnes { 356dfbe8321SBarry Smith PetscErrorCode ierr; 357ace3abfcSBarry Smith PetscBool iascii; 358a935fc98SLois Curfman McInnes 3593a40ed3dSBarry Smith PetscFunctionBegin; 360251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 36132077d6dSBarry Smith if (iascii) { 36250f6de3fSJed Brown } 36350f6de3fSJed Brown PetscFunctionReturn(0); 36450f6de3fSJed Brown } 36550f6de3fSJed Brown 36604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 36704d965bbSLois Curfman McInnes /* 36804d7464bSBarry Smith SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method. 36904d965bbSLois Curfman McInnes 37004d965bbSLois Curfman McInnes Input Parameter: 37104d965bbSLois Curfman McInnes . snes - the SNES context 37204d965bbSLois Curfman McInnes 37304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 37404d965bbSLois Curfman McInnes */ 3754a2ae208SSatish Balay #undef __FUNCT__ 37604d7464bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NEWTONLS" 3778c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_NEWTONLS(PetscOptions *PetscOptionsObject,SNES snes) 3785e42470aSBarry Smith { 379dfbe8321SBarry Smith PetscErrorCode ierr; 380f1c6b773SPeter Brune SNESLineSearch linesearch; 3815e42470aSBarry Smith 3823a40ed3dSBarry Smith PetscFunctionBegin; 3839e764e56SPeter Brune if (!snes->linesearch) { 3847601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3851a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);CHKERRQ(ierr); 3869e764e56SPeter Brune } 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 3885e42470aSBarry Smith } 3894827ddcaSBarry Smith 39004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 391ebe3b25bSBarry Smith /*MC 39204d7464bSBarry Smith SNESNEWTONLS - Newton based nonlinear solver that uses a line search 39304d965bbSLois Curfman McInnes 394ebe3b25bSBarry Smith Options Database: 3955a9b6599SPeter Brune + -snes_linesearch_type <bt> - bt,basic. Select line search type 3965a9b6599SPeter Brune . -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt 3975a9b6599SPeter Brune . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch 39882d26c24SPeter Brune . -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient 39982d26c24SPeter 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) 4005a9b6599SPeter Brune . -snes_linesearch_minlambda <minlambda> - Sets the minimum lambda the line search will tolerate 40182d26c24SPeter Brune . -snes_linesearch_monitor - print information about progress of line searches 40282d26c24SPeter Brune - -snes_linesearch_damping - damping factor used for basic line search 403acbee50cSBarry Smith 404ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 40504d965bbSLois Curfman McInnes 406ee3001cbSBarry Smith Level: beginner 407ee3001cbSBarry Smith 40804d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder() 4095a9b6599SPeter Brune SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms() 410ebe3b25bSBarry Smith 411ebe3b25bSBarry Smith M*/ 4124a2ae208SSatish Balay #undef __FUNCT__ 41304d7464bSBarry Smith #define __FUNCT__ "SNESCreate_NEWTONLS" 4148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) 4155e42470aSBarry Smith { 416dfbe8321SBarry Smith PetscErrorCode ierr; 41704d7464bSBarry Smith SNES_NEWTONLS *neP; 4185e42470aSBarry Smith 4193a40ed3dSBarry Smith PetscFunctionBegin; 42004d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONLS; 42104d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONLS; 42204d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONLS; 42304d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS; 42404d7464bSBarry Smith snes->ops->view = SNESView_NEWTONLS; 42504d7464bSBarry Smith snes->ops->reset = SNESReset_NEWTONLS; 4265e42470aSBarry Smith 4276c67d002SPeter Brune snes->pcside = PC_RIGHT; 42842f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 4296b6c8463SPeter Brune snes->usespc = PETSC_TRUE; 430b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 4315e42470aSBarry Smith snes->data = (void*)neP; 4323a40ed3dSBarry Smith PetscFunctionReturn(0); 4335e42470aSBarry Smith } 434