163dd3a1aSKris Buschelman #define PETSCSNES_DLL 25e76c431SBarry Smith 370f55243SBarry Smith #include "src/snes/impls/ls/ls.h" 45e42470aSBarry Smith 58a5d9ee4SBarry Smith /* 68a5d9ee4SBarry Smith Checks if J^T F = 0 which implies we've found a local minimum of the function, 78a5d9ee4SBarry Smith but not a zero. In the case when one cannot compute J^T F we use the fact that 836669109SBarry Smith 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 936669109SBarry Smith for this trick. 108a5d9ee4SBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private" 13dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 148a5d9ee4SBarry Smith { 158a5d9ee4SBarry Smith PetscReal a1; 16dfbe8321SBarry Smith PetscErrorCode ierr; 1736669109SBarry Smith PetscTruth hastranspose; 188a5d9ee4SBarry Smith 198a5d9ee4SBarry Smith PetscFunctionBegin; 208a5d9ee4SBarry Smith *ismin = PETSC_FALSE; 2136669109SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 2236669109SBarry Smith if (hastranspose) { 238a5d9ee4SBarry Smith /* Compute || J^T F|| */ 248a5d9ee4SBarry Smith ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 258a5d9ee4SBarry Smith ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 2663ba0a88SBarry Smith ierr = PetscLogInfo((0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm));CHKERRQ(ierr); 278a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2836669109SBarry Smith } else { 2936669109SBarry Smith Vec work; 30ea709b57SSatish Balay PetscScalar result; 3136669109SBarry Smith PetscReal wnorm; 3236669109SBarry Smith 3336669109SBarry Smith ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr); 3436669109SBarry Smith ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 3536669109SBarry Smith ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 3636669109SBarry Smith ierr = MatMult(A,W,work);CHKERRQ(ierr); 3736669109SBarry Smith ierr = VecDot(F,work,&result);CHKERRQ(ierr); 3836669109SBarry Smith ierr = VecDestroy(work);CHKERRQ(ierr); 3936669109SBarry Smith a1 = PetscAbsScalar(result)/(fnorm*wnorm); 4063ba0a88SBarry Smith ierr = PetscLogInfo((0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1));CHKERRQ(ierr); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4674637425SBarry Smith /* 475ed2d596SBarry Smith Checks if J^T(F - J*X) = 0 4874637425SBarry Smith */ 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private" 51dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 5274637425SBarry Smith { 5374637425SBarry Smith PetscReal a1,a2; 54dfbe8321SBarry Smith PetscErrorCode ierr; 5574637425SBarry Smith PetscTruth hastranspose; 56ea709b57SSatish Balay PetscScalar mone = -1.0; 5774637425SBarry Smith 5874637425SBarry Smith PetscFunctionBegin; 5974637425SBarry Smith ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 6074637425SBarry Smith if (hastranspose) { 6174637425SBarry Smith ierr = MatMult(A,X,W1);CHKERRQ(ierr); 6274637425SBarry Smith ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr); 6374637425SBarry Smith 6474637425SBarry Smith /* Compute || J^T W|| */ 6574637425SBarry Smith ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 6674637425SBarry Smith ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 6774637425SBarry Smith ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 689994e62eSSatish Balay if (a1 != 0) { 6963ba0a88SBarry Smith ierr = PetscLogInfo((0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1));CHKERRQ(ierr); 7085471664SBarry Smith } 7174637425SBarry Smith } 7274637425SBarry Smith PetscFunctionReturn(0); 7374637425SBarry Smith } 7474637425SBarry Smith 7504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 7604d965bbSLois Curfman McInnes 7704d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 7894b7f48cSBarry Smith for solving a system of nonlinear equations, using the KSP, Vec, 7904d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 8004d965bbSLois Curfman McInnes respectively. 8104d965bbSLois Curfman McInnes 82fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 8304d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 8404d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 85fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 8604d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 8704d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 884b27c08aSLois Curfman McInnes we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 894b27c08aSLois Curfman McInnes systems of nonlinear equations with a line search (LS) method. 9004d965bbSLois Curfman McInnes These routines are actually called via the common user interface 9104d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 9204d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 9304d965bbSLois Curfman McInnes for all nonlinear solvers. 9404d965bbSLois Curfman McInnes 9504d965bbSLois Curfman McInnes Another key routine is: 9604d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 9704d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 9804d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 9904d965bbSLois Curfman McInnes SNESSolve() if necessary. 10004d965bbSLois Curfman McInnes 10104d965bbSLois Curfman McInnes Additional basic routines are: 10204d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 10304d965bbSLois Curfman McInnes have actually been used. 10404d965bbSLois Curfman McInnes These are called by application codes via the interface routines 105186905e3SBarry Smith SNESView(). 10604d965bbSLois Curfman McInnes 10704d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 10804d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 10904d965bbSLois Curfman McInnes above description applies to these categories also. 11004d965bbSLois Curfman McInnes 11104d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 1125e42470aSBarry Smith /* 1134b27c08aSLois Curfman McInnes SNESSolve_LS - Solves a nonlinear system with a truncated Newton 11404d965bbSLois Curfman McInnes method with a line search. 1155e76c431SBarry Smith 11604d965bbSLois Curfman McInnes Input Parameters: 11704d965bbSLois Curfman McInnes . snes - the SNES context 1185e76c431SBarry Smith 11904d965bbSLois Curfman McInnes Output Parameter: 120c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 12104d965bbSLois Curfman McInnes 12204d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 1235e76c431SBarry Smith 1245e76c431SBarry Smith Notes: 1255e76c431SBarry Smith This implements essentially a truncated Newton method with a 1265e76c431SBarry Smith line search. By default a cubic backtracking line search 1275e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1285e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 129393d2d9aSLois Curfman McInnes and Schnabel. 1305e76c431SBarry Smith */ 1314a2ae208SSatish Balay #undef __FUNCT__ 1324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS" 133dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes) 1345e76c431SBarry Smith { 1354b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 1366849ba73SBarry Smith PetscErrorCode ierr; 137a7cc72afSBarry Smith PetscInt maxits,i,lits; 138a7cc72afSBarry Smith PetscTruth lssucceed; 139112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 140329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1415e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 142c293cc10SBarry Smith KSP ksp; 1435e76c431SBarry Smith 1443a40ed3dSBarry Smith PetscFunctionBegin; 14594b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 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 */ 1515e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1525e42470aSBarry Smith G = snes->work[1]; 1535e42470aSBarry Smith W = snes->work[2]; 1545e76c431SBarry Smith 1554c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1564c49b128SBarry Smith snes->iter = 0; 1574c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 158*19717074SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 159cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 16043e71028SBarry Smith if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 1610f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1625e42470aSBarry Smith snes->norm = fnorm; 1630f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 16484c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 16594a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1665e76c431SBarry Smith 16770441072SBarry Smith if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 16894a9d846SBarry Smith 169d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 170d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 171d034289bSBarry Smith 1725e76c431SBarry Smith for (i=0; i<maxits; i++) { 1735e76c431SBarry Smith 174329e7e40SMatthew Knepley /* Call general purpose update function */ 175abc0a331SBarry Smith if (snes->update) { 176329e7e40SMatthew Knepley ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr); 177de8cb200SMatthew Knepley } 178329e7e40SMatthew Knepley 179ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1805334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 18194b7f48cSBarry Smith ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 18223ce1328SBarry Smith ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 183c293cc10SBarry Smith ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr); 18474637425SBarry Smith 185b0a32e0cSBarry Smith if (PetscLogPrintInfo){ 18674637425SBarry Smith ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 18785471664SBarry Smith } 18874637425SBarry Smith 18990cbd584SBarry Smith /* should check what happened to the linear solve? */ 1903505fcc1SBarry Smith snes->linear_its += lits; 19163ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr); 192ea4d3ed3SLois Curfman McInnes 193ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 194ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 195ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 196ea4d3ed3SLois Curfman McInnes */ 19781b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 198a7cc72afSBarry Smith ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 19963ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed));CHKERRQ(ierr); 20043e71028SBarry Smith if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 201a3b891d8SBarry Smith 202a3b891d8SBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 203a3b891d8SBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 204a3b891d8SBarry Smith fnorm = gnorm; 205a3b891d8SBarry Smith 2065ed2d596SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2075ed2d596SBarry Smith snes->iter = i+1; 2085ed2d596SBarry Smith snes->norm = fnorm; 2095ed2d596SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2105ed2d596SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 2115ed2d596SBarry Smith SNESMonitor(snes,i+1,fnorm); 2125ed2d596SBarry Smith 213a7cc72afSBarry Smith if (!lssucceed) { 2148a5d9ee4SBarry Smith PetscTruth ismin; 21550ffb88aSMatthew Knepley 21650ffb88aSMatthew Knepley if (++snes->numFailures >= snes->maxFailures) { 2173505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 2188a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 2198a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 2203505fcc1SBarry Smith break; 2213505fcc1SBarry Smith } 22250ffb88aSMatthew Knepley } 2235e76c431SBarry Smith 2245e76c431SBarry Smith /* Test for convergence */ 22529e0b56fSBarry Smith if (snes->converged) { 22629e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 2273505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2283505fcc1SBarry Smith if (snes->reason) { 22916c95cacSBarry Smith break; 23016c95cacSBarry Smith } 23116c95cacSBarry Smith } 23229e0b56fSBarry Smith } 23339e2f89bSBarry Smith if (X != snes->vec_sol) { 234393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 235e15f7bb5SBarry Smith } 236e15f7bb5SBarry Smith if (F != snes->vec_func) { 237e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 238e15f7bb5SBarry Smith } 23939e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 24039e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 24152392280SLois Curfman McInnes if (i == maxits) { 24263ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr); 2433505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 24452392280SLois Curfman McInnes } 2453a40ed3dSBarry Smith PetscFunctionReturn(0); 2465e76c431SBarry Smith } 24704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 24804d965bbSLois Curfman McInnes /* 2494b27c08aSLois Curfman McInnes SNESSetUp_LS - Sets up the internal data structures for the later use 2504b27c08aSLois Curfman McInnes of the SNESLS nonlinear solver. 25104d965bbSLois Curfman McInnes 25204d965bbSLois Curfman McInnes Input Parameter: 25304d965bbSLois Curfman McInnes . snes - the SNES context 25404d965bbSLois Curfman McInnes . x - the solution vector 25504d965bbSLois Curfman McInnes 25604d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 25704d965bbSLois Curfman McInnes 25804d965bbSLois Curfman McInnes Notes: 2594b27c08aSLois Curfman McInnes For basic use of the SNES solvers, the user need not explicitly call 26004d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 26104d965bbSLois Curfman McInnes the call to SNESSolve(). 26204d965bbSLois Curfman McInnes */ 2634a2ae208SSatish Balay #undef __FUNCT__ 2644b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS" 265dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes) 2665e76c431SBarry Smith { 267dfbe8321SBarry Smith PetscErrorCode ierr; 2683a40ed3dSBarry Smith 2693a40ed3dSBarry Smith PetscFunctionBegin; 27081b6cf68SLois Curfman McInnes snes->nwork = 4; 271d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 272efee365bSSatish Balay ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 27381b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2743a40ed3dSBarry Smith PetscFunctionReturn(0); 2755e76c431SBarry Smith } 27604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 27704d965bbSLois Curfman McInnes /* 2784b27c08aSLois Curfman McInnes SNESDestroy_LS - Destroys the private SNES_LS context that was created 2794b27c08aSLois Curfman McInnes with SNESCreate_LS(). 28004d965bbSLois Curfman McInnes 28104d965bbSLois Curfman McInnes Input Parameter: 28204d965bbSLois Curfman McInnes . snes - the SNES context 28304d965bbSLois Curfman McInnes 284de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 28504d965bbSLois Curfman McInnes */ 2864a2ae208SSatish Balay #undef __FUNCT__ 2874b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS" 288dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes) 2895e76c431SBarry Smith { 290dfbe8321SBarry Smith PetscErrorCode ierr; 2913a40ed3dSBarry Smith 2923a40ed3dSBarry Smith PetscFunctionBegin; 2935baf8537SBarry Smith if (snes->nwork) { 2944b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 29521c89e3eSBarry Smith } 2965c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2973a40ed3dSBarry Smith PetscFunctionReturn(0); 2985e76c431SBarry Smith } 29904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3004a2ae208SSatish Balay #undef __FUNCT__ 3014a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch" 30282bf6240SBarry Smith 3034b828684SBarry Smith /*@C 3045e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 3055e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 3065e76c431SBarry Smith to serve as a template and is not recommended for general use. 3075e76c431SBarry Smith 308c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 309c7afd0dbSLois Curfman McInnes 3105e76c431SBarry Smith Input Parameters: 311c7afd0dbSLois Curfman McInnes + snes - nonlinear context 312af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 3135e76c431SBarry Smith . x - current iterate 3145e76c431SBarry Smith . f - residual evaluated at x 3155e76c431SBarry Smith . y - search direction (contains new iterate on output) 3165e76c431SBarry Smith . w - work vector 317c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 3185e76c431SBarry Smith 319c4a48953SLois Curfman McInnes Output Parameters: 320c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 3215e76c431SBarry Smith . y - new iterate (contains search direction on input) 3225e76c431SBarry Smith . gnorm - 2-norm of g 3235e76c431SBarry Smith . ynorm - 2-norm of search length 324a7cc72afSBarry Smith - flag - PETSC_TRUE on success, PETSC_FALSE on failure 325fee21e36SBarry Smith 326c4a48953SLois Curfman McInnes Options Database Key: 3274b27c08aSLois Curfman McInnes . -snes_ls basic - Activates SNESNoLineSearch() 328c4a48953SLois Curfman McInnes 32936851e7fSLois Curfman McInnes Level: advanced 33036851e7fSLois Curfman McInnes 33128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 33228ae5a21SLois Curfman McInnes 333f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 334af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 3355e76c431SBarry Smith @*/ 33663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 3375e76c431SBarry Smith { 338dfbe8321SBarry Smith PetscErrorCode ierr; 339ea709b57SSatish Balay PetscScalar mone = -1.0; 3404b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 3418f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 3425334005bSBarry Smith 3433a40ed3dSBarry Smith PetscFunctionBegin; 344a7cc72afSBarry Smith *flag = PETSC_TRUE; 345d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 346a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 347ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3488f99978dSLois Curfman McInnes if (neP->CheckStep) { 3498f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3508f99978dSLois Curfman McInnes } 351*19717074SBarry Smith ierr = SNESComputeFunction(snes,y,g); 352*19717074SBarry Smith if (ierr) { 353*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 354*19717074SBarry Smith CHKERRQ(ierr); 355*19717074SBarry Smith } 356a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 35743e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 358d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 3593a40ed3dSBarry Smith PetscFunctionReturn(0); 3605e76c431SBarry Smith } 36104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3624a2ae208SSatish Balay #undef __FUNCT__ 3634a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms" 36482bf6240SBarry Smith 36529e0b56fSBarry Smith /*@C 36629e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 36729e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 36829e0b56fSBarry Smith even compute the norm of the function or search direction; this 36929e0b56fSBarry Smith is intended only when you know the full step is fine and are 37029e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 371c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 372c7afd0dbSLois Curfman McInnes 373c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 37429e0b56fSBarry Smith 37529e0b56fSBarry Smith Input Parameters: 376c7afd0dbSLois Curfman McInnes + snes - nonlinear context 377af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 37829e0b56fSBarry Smith . x - current iterate 37929e0b56fSBarry Smith . f - residual evaluated at x 38029e0b56fSBarry Smith . y - search direction (contains new iterate on output) 38129e0b56fSBarry Smith . w - work vector 382c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 38329e0b56fSBarry Smith 38429e0b56fSBarry Smith Output Parameters: 385c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 38629e0b56fSBarry Smith . gnorm - not changed 38729e0b56fSBarry Smith . ynorm - not changed 388a7cc72afSBarry Smith - flag - set to PETSC_TRUE indicating a successful line search 389fee21e36SBarry Smith 39029e0b56fSBarry Smith Options Database Key: 3914b27c08aSLois Curfman McInnes . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 39229e0b56fSBarry Smith 3938cbba510SBarry Smith Notes: 394ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 395ea56f5baSLois Curfman McInnes either the options 396ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 397ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 398329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 399329f5518SBarry Smith otherwise, the SNES solver will generate an error. 400329f5518SBarry Smith 401329f5518SBarry Smith During the final iteration this will not evaluate the function at 402329f5518SBarry Smith the solution point. This is to save a function evaluation while 403329f5518SBarry Smith using pseudo-timestepping. 4048cbba510SBarry Smith 405ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 406ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 407ea56f5baSLois Curfman McInnes correct, since they are not computed. 408ea56f5baSLois Curfman McInnes 409ea56f5baSLois Curfman McInnes Level: advanced 4108cbba510SBarry Smith 41129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 41229e0b56fSBarry Smith 41329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 41429e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 41529e0b56fSBarry Smith @*/ 41663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 41729e0b56fSBarry Smith { 418dfbe8321SBarry Smith PetscErrorCode ierr; 419ea709b57SSatish Balay PetscScalar mone = -1.0; 4204b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 4218f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 42229e0b56fSBarry Smith 4233a40ed3dSBarry Smith PetscFunctionBegin; 424a7cc72afSBarry Smith *flag = PETSC_TRUE; 425d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 42629e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 4278f99978dSLois Curfman McInnes if (neP->CheckStep) { 4288f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 4298f99978dSLois Curfman McInnes } 430329f5518SBarry Smith 431329f5518SBarry Smith /* don't evaluate function the last time through */ 432329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 433*19717074SBarry Smith ierr = SNESComputeFunction(snes,y,g); 434*19717074SBarry Smith if (ierr) { 435*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 436*19717074SBarry Smith CHKERRQ(ierr); 437*19717074SBarry Smith } 438329f5518SBarry Smith } 439d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 4403a40ed3dSBarry Smith PetscFunctionReturn(0); 44129e0b56fSBarry Smith } 44204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 4434a2ae208SSatish Balay #undef __FUNCT__ 4444a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch" 4454b828684SBarry Smith /*@C 446f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 4475e76c431SBarry Smith 448c7afd0dbSLois Curfman McInnes Collective on SNES 449c7afd0dbSLois Curfman McInnes 4505e76c431SBarry Smith Input Parameters: 451c7afd0dbSLois Curfman McInnes + snes - nonlinear context 452af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4535e76c431SBarry Smith . x - current iterate 4545e76c431SBarry Smith . f - residual evaluated at x 4555e76c431SBarry Smith . y - search direction (contains new iterate on output) 4565e76c431SBarry Smith . w - work vector 457c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4585e76c431SBarry Smith 459393d2d9aSLois Curfman McInnes Output Parameters: 460c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4615e76c431SBarry Smith . y - new iterate (contains search direction on input) 4625e76c431SBarry Smith . gnorm - 2-norm of g 4635e76c431SBarry Smith . ynorm - 2-norm of search length 464a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 465fee21e36SBarry Smith 466c4a48953SLois Curfman McInnes Options Database Key: 4674b27c08aSLois Curfman McInnes $ -snes_ls cubic - Activates SNESCubicLineSearch() 468c4a48953SLois Curfman McInnes 4695e76c431SBarry Smith Notes: 4705e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4715e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4725e76c431SBarry Smith 47336851e7fSLois Curfman McInnes Level: advanced 47436851e7fSLois Curfman McInnes 47528ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 47628ae5a21SLois Curfman McInnes 477af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4785e76c431SBarry Smith @*/ 47963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 4805e76c431SBarry Smith { 481406556e6SLois Curfman McInnes /* 482406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 483406556e6SLois Curfman McInnes minimization problem: 484406556e6SLois Curfman McInnes min z(x): R^n -> R, 485406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 486406556e6SLois Curfman McInnes */ 487406556e6SLois Curfman McInnes 4885ca10a37SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 489329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 490aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 491ea709b57SSatish Balay PetscScalar cinitslope,clambda; 4926b5873e3SBarry Smith #endif 493dfbe8321SBarry Smith PetscErrorCode ierr; 494a7cc72afSBarry Smith PetscInt count; 4954b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 496ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 4978f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4985e76c431SBarry Smith 4993a40ed3dSBarry Smith PetscFunctionBegin; 500d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 501a7cc72afSBarry Smith *flag = PETSC_TRUE; 5025e76c431SBarry Smith alpha = neP->alpha; 5035e76c431SBarry Smith maxstep = neP->maxstep; 5045e76c431SBarry Smith steptol = neP->steptol; 5055e76c431SBarry Smith 506cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 50763b9fa5eSBarry Smith if (*ynorm == 0.0) { 50863ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr); 509a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 510a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 511a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 512a7cc72afSBarry Smith *flag = PETSC_FALSE; 513ad922ac9SBarry Smith goto theend1; 51494a9d846SBarry Smith } 5155e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 5165e42470aSBarry Smith scale = maxstep/(*ynorm); 517aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 51863ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr); 5196b5873e3SBarry Smith #else 52063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr); 5216b5873e3SBarry Smith #endif 522393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 5235e76c431SBarry Smith *ynorm = maxstep; 5245e76c431SBarry Smith } 525ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 5265ca10a37SBarry Smith minlambda = steptol/rellength; 527a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 528aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 529a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 530329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 5315e42470aSBarry Smith #else 532a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 5335e42470aSBarry Smith #endif 5345e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 5355e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 5365e76c431SBarry Smith 537393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 5385334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 53943e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 54063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr); 54143e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 54243e71028SBarry Smith *flag = PETSC_FALSE; 54343e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 54443e71028SBarry Smith goto theend1; 54543e71028SBarry Smith } 546*19717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 547*19717074SBarry Smith if (ierr) { 548*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 549*19717074SBarry Smith CHKERRQ(ierr); 550*19717074SBarry Smith } 551313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 55243e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5535d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 554393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 55563ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Using full step\n"));CHKERRQ(ierr); 556ad922ac9SBarry Smith goto theend1; 5575e76c431SBarry Smith } 5585e76c431SBarry Smith 5595e76c431SBarry Smith /* Fit points with quadratic */ 560313b5538SBarry Smith lambda = 1.0; 561a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 5625e76c431SBarry Smith lambdaprev = lambda; 56343e71028SBarry Smith printf("tmp %g ddsd %g %g %g %g\n",lambdatemp,((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope),*gnorm,fnorm,initslope); 5645e76c431SBarry Smith gnormprev = *gnorm; 565329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 566ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 567ddd12767SBarry Smith else lambda = lambdatemp; 568393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 569ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 570aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 571ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5725e42470aSBarry Smith #else 573ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5745e42470aSBarry Smith #endif 57543e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 57663ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr); 57743e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 57843e71028SBarry Smith *flag = PETSC_FALSE; 57943e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 58043e71028SBarry Smith goto theend1; 58143e71028SBarry Smith } 582*19717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 583*19717074SBarry Smith if (ierr) { 584*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 585*19717074SBarry Smith CHKERRQ(ierr); 586*19717074SBarry Smith } 587cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 58843e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 5895ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 590393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 59163ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 592ad922ac9SBarry Smith goto theend1; 5935e76c431SBarry Smith } 5945e76c431SBarry Smith 5955e76c431SBarry Smith /* Fit points with cubic */ 5965e76c431SBarry Smith count = 1; 5978229bfc2SMatthew Knepley while (count < 10000) { 5985e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 59963ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr); 60063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 601f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 60243e71028SBarry Smith *flag = PETSC_FALSE; 60343e71028SBarry Smith break; 6045e76c431SBarry Smith } 60543e71028SBarry Smith printf("lamd %g %g\n",lambda,lambdaprev); 6065d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 6075d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 608ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6092b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 6105e76c431SBarry Smith d = b*b - 3*a*initslope; 6115e76c431SBarry Smith if (d < 0.0) d = 0.0; 6125e76c431SBarry Smith if (a == 0.0) { 6135e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 6145e76c431SBarry Smith } else { 6155e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 6165e76c431SBarry Smith } 6175e76c431SBarry Smith lambdaprev = lambda; 6185e76c431SBarry Smith gnormprev = *gnorm; 619329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 620329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 6215e76c431SBarry Smith else lambda = lambdatemp; 622393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 623ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 624aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 625ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 626393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6275e42470aSBarry Smith #else 628ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 6295e42470aSBarry Smith #endif 63043e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 63163ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr); 63263ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 63343e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 634ed98166eSMatthew Knepley *flag = PETSC_FALSE; 63543e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 636ed98166eSMatthew Knepley break; 637ed98166eSMatthew Knepley } 638*19717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 639*19717074SBarry Smith if (ierr) { 640*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 641*19717074SBarry Smith CHKERRQ(ierr); 642*19717074SBarry Smith } 643cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 64443e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 6455ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 646393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 64763ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 64843e71028SBarry Smith break; 6492b022350SLois Curfman McInnes } else { 65063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 6515e76c431SBarry Smith } 6525e76c431SBarry Smith count++; 6535e76c431SBarry Smith } 6548229bfc2SMatthew Knepley if (count >= 10000) { 655cd7625f8SBarry Smith SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 6568229bfc2SMatthew Knepley } 657ad922ac9SBarry Smith theend1: 6588f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 65943e71028SBarry Smith if (neP->CheckStep && *flag) { 6608f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 661abc0a331SBarry Smith if (change_y) { /* recompute the function if the step has changed */ 662*19717074SBarry Smith ierr = SNESComputeFunction(snes,y,g); 663*19717074SBarry Smith if (ierr) { 664*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 665*19717074SBarry Smith CHKERRQ(ierr); 666*19717074SBarry Smith } 6678f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 66843e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 66943e71028SBarry Smith ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 6708f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 6718f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 6728f99978dSLois Curfman McInnes } 6738f99978dSLois Curfman McInnes } 674d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 6753a40ed3dSBarry Smith PetscFunctionReturn(0); 6765e76c431SBarry Smith } 67704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 6784a2ae208SSatish Balay #undef __FUNCT__ 6794a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch" 6804b828684SBarry Smith /*@C 681f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 6825e76c431SBarry Smith 683c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 684c7afd0dbSLois Curfman McInnes 6855e42470aSBarry Smith Input Parameters: 686c7afd0dbSLois Curfman McInnes + snes - the SNES context 687af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 6885e42470aSBarry Smith . x - current iterate 6895e42470aSBarry Smith . f - residual evaluated at x 6905e42470aSBarry Smith . y - search direction (contains new iterate on output) 6915e42470aSBarry Smith . w - work vector 692c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 6935e42470aSBarry Smith 694c4a48953SLois Curfman McInnes Output Parameters: 695c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 6965e42470aSBarry Smith . y - new iterate (contains search direction on input) 6975e42470aSBarry Smith . gnorm - 2-norm of g 6985e42470aSBarry Smith . ynorm - 2-norm of search length 699a7cc72afSBarry Smith - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 700fee21e36SBarry Smith 701c4a48953SLois Curfman McInnes Options Database Key: 7024b27c08aSLois Curfman McInnes . -snes_ls quadratic - Activates SNESQuadraticLineSearch() 703c4a48953SLois Curfman McInnes 7045e42470aSBarry Smith Notes: 7054b27c08aSLois Curfman McInnes Use SNESSetLineSearch() to set this routine within the SNESLS method. 7065e42470aSBarry Smith 70736851e7fSLois Curfman McInnes Level: advanced 70836851e7fSLois Curfman McInnes 709f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 710f59ffdedSLois Curfman McInnes 711af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 7125e42470aSBarry Smith @*/ 71363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 7145e76c431SBarry Smith { 715406556e6SLois Curfman McInnes /* 716406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 717406556e6SLois Curfman McInnes minimization problem: 718406556e6SLois Curfman McInnes min z(x): R^n -> R, 719406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 720406556e6SLois Curfman McInnes */ 7215ca10a37SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 722aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 723ea709b57SSatish Balay PetscScalar cinitslope,clambda; 7246b5873e3SBarry Smith #endif 725dfbe8321SBarry Smith PetscErrorCode ierr; 726a7cc72afSBarry Smith PetscInt count; 7274b27c08aSLois Curfman McInnes SNES_LS *neP = (SNES_LS*)snes->data; 728ea709b57SSatish Balay PetscScalar mone = -1.0,scale; 7298f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 7305e76c431SBarry Smith 7313a40ed3dSBarry Smith PetscFunctionBegin; 732d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 733a7cc72afSBarry Smith *flag = PETSC_TRUE; 7345e42470aSBarry Smith alpha = neP->alpha; 7355e42470aSBarry Smith maxstep = neP->maxstep; 7365e42470aSBarry Smith steptol = neP->steptol; 7375e76c431SBarry Smith 7383505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 73963b9fa5eSBarry Smith if (*ynorm == 0.0) { 74063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr); 741b37302c6SLois Curfman McInnes *gnorm = fnorm; 742b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 743b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 744a7cc72afSBarry Smith *flag = PETSC_FALSE; 745ad922ac9SBarry Smith goto theend2; 74694a9d846SBarry Smith } 7475e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 7485e42470aSBarry Smith scale = maxstep/(*ynorm); 749393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 7505e42470aSBarry Smith *ynorm = maxstep; 7515e76c431SBarry Smith } 752ba6a83e5SMatthew Knepley ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 7535ca10a37SBarry Smith minlambda = steptol/rellength; 754a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 755aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 756a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 757329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 7585e42470aSBarry Smith #else 759a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 7605e42470aSBarry Smith #endif 7615e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 7625e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 7635e42470aSBarry Smith 764393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 7655334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 76643e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 76763ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr); 76843e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 76943e71028SBarry Smith *flag = PETSC_FALSE; 77043e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 77143e71028SBarry Smith goto theend2; 77243e71028SBarry Smith } 773*19717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 774*19717074SBarry Smith if (ierr) { 775*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 776*19717074SBarry Smith CHKERRQ(ierr); 777*19717074SBarry Smith } 778cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 77943e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 7805d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 781393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 78263ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Using full step\n"));CHKERRQ(ierr); 783ad922ac9SBarry Smith goto theend2; 7845e42470aSBarry Smith } 7855e42470aSBarry Smith 7865e42470aSBarry Smith /* Fit points with quadratic */ 787313b5538SBarry Smith lambda = 1.0; 7885e42470aSBarry Smith count = 1; 7895ca10a37SBarry Smith while (PETSC_TRUE) { 7905e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 79163ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr); 79263ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 793f168f371SBarry Smith ierr = VecCopy(x,y);CHKERRQ(ierr); 79443e71028SBarry Smith *flag = PETSC_FALSE; 79543e71028SBarry Smith break; 7965e42470aSBarry Smith } 797a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 798329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 799329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 800329f5518SBarry Smith else lambda = lambdatemp; 801393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 8023505fcc1SBarry Smith lambdaneg = -lambda; 803aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 8043505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 8055e42470aSBarry Smith #else 8063505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 8075e42470aSBarry Smith #endif 80843e71028SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 80963ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr); 81063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 81143e71028SBarry Smith ierr = VecCopy(w,y);CHKERRQ(ierr); 812ed98166eSMatthew Knepley *flag = PETSC_FALSE; 81343e71028SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 814ed98166eSMatthew Knepley break; 815ed98166eSMatthew Knepley } 816*19717074SBarry Smith ierr = SNESComputeFunction(snes,w,g); 817*19717074SBarry Smith if (ierr) { 818*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 819*19717074SBarry Smith CHKERRQ(ierr); 820*19717074SBarry Smith } 821cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 82243e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 8235ed2d596SBarry Smith if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 824393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 82563ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr); 82606259719SBarry Smith break; 8275e42470aSBarry Smith } 8285e42470aSBarry Smith count++; 8295e42470aSBarry Smith } 830ad922ac9SBarry Smith theend2: 8318f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 8328f99978dSLois Curfman McInnes if (neP->CheckStep) { 8338f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 834abc0a331SBarry Smith if (change_y) { /* recompute the function if the step has changed */ 835*19717074SBarry Smith ierr = SNESComputeFunction(snes,y,g); 836*19717074SBarry Smith if (ierr) { 837*19717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 838*19717074SBarry Smith CHKERRQ(ierr); 839*19717074SBarry Smith } 8408f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 84143e71028SBarry Smith if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 84243e71028SBarry Smith ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 8438f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 8448f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 8458f99978dSLois Curfman McInnes } 8468f99978dSLois Curfman McInnes } 847d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 8483a40ed3dSBarry Smith PetscFunctionReturn(0); 8495e76c431SBarry Smith } 8502343ba6eSBarry Smith 85104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 8524a2ae208SSatish Balay #undef __FUNCT__ 8534a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch" 854c9e6a524SLois Curfman McInnes /*@C 85577c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 8564b27c08aSLois Curfman McInnes by the method SNESLS. 8575e76c431SBarry Smith 8585e76c431SBarry Smith Input Parameters: 859c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 860af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 861c7afd0dbSLois Curfman McInnes - func - pointer to int function 8625e76c431SBarry Smith 863fee21e36SBarry Smith Collective on SNES 864fee21e36SBarry Smith 865c4a48953SLois Curfman McInnes Available Routines: 866c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 867f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 868f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 869af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 8705e76c431SBarry Smith 871c4a48953SLois Curfman McInnes Options Database Keys: 8724b27c08aSLois Curfman McInnes + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 8734b27c08aSLois Curfman McInnes . -snes_ls_alpha <alpha> - Sets alpha 8744b27c08aSLois Curfman McInnes . -snes_ls_maxstep <max> - Sets maxstep 8754b27c08aSLois Curfman McInnes - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 8763304466cSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 8773304466cSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 878c4a48953SLois Curfman McInnes 8795e76c431SBarry Smith Calling sequence of func: 880bd208895SLois Curfman McInnes .vb 881af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 882a7cc72afSBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 883bd208895SLois Curfman McInnes .ve 8845e76c431SBarry Smith 8855e76c431SBarry Smith Input parameters for func: 886c7afd0dbSLois Curfman McInnes + snes - nonlinear context 887af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 8885e76c431SBarry Smith . x - current iterate 8895e76c431SBarry Smith . f - residual evaluated at x 8905e76c431SBarry Smith . y - search direction (contains new iterate on output) 8915e76c431SBarry Smith . w - work vector 892c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 8935e76c431SBarry Smith 8945e76c431SBarry Smith Output parameters for func: 895c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 8965e76c431SBarry Smith . y - new iterate (contains search direction on input) 8975e76c431SBarry Smith . gnorm - 2-norm of g 8985e76c431SBarry Smith . ynorm - 2-norm of search length 899a7cc72afSBarry Smith - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 900f59ffdedSLois Curfman McInnes 90136851e7fSLois Curfman McInnes Level: advanced 90236851e7fSLois Curfman McInnes 903f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 904f59ffdedSLois Curfman McInnes 9055d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 9065d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 9075e76c431SBarry Smith @*/ 90863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 9095e76c431SBarry Smith { 910a7cc72afSBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 91182bf6240SBarry Smith 9123a40ed3dSBarry Smith PetscFunctionBegin; 913c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr); 91482bf6240SBarry Smith if (f) { 915af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 91682bf6240SBarry Smith } 9173a40ed3dSBarry Smith PetscFunctionReturn(0); 9185e76c431SBarry Smith } 9198e019c35SBarry Smith 920a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 92104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 922fb2e594dSBarry Smith EXTERN_C_BEGIN 9234a2ae208SSatish Balay #undef __FUNCT__ 9244a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS" 92563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx) 92682bf6240SBarry Smith { 92782bf6240SBarry Smith PetscFunctionBegin; 9284b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->LineSearch = func; 9294b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->lsP = lsctx; 93082bf6240SBarry Smith PetscFunctionReturn(0); 93182bf6240SBarry Smith } 932fb2e594dSBarry Smith EXTERN_C_END 93304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9344a2ae208SSatish Balay #undef __FUNCT__ 9354a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck" 936c8dd0c0dSLois Curfman McInnes /*@C 937530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 9384b27c08aSLois Curfman McInnes by the line search routine in the Newton-based method SNESLS. 939c8dd0c0dSLois Curfman McInnes 940c8dd0c0dSLois Curfman McInnes Input Parameters: 941c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 942c8dd0c0dSLois Curfman McInnes . func - pointer to int function 943c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 944c8dd0c0dSLois Curfman McInnes 945c8dd0c0dSLois Curfman McInnes Collective on SNES 946c8dd0c0dSLois Curfman McInnes 947c8dd0c0dSLois Curfman McInnes Calling sequence of func: 948c8dd0c0dSLois Curfman McInnes .vb 949b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 950c8dd0c0dSLois Curfman McInnes .ve 951b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 952b1ae27eaSLois Curfman McInnes on failure. 953c8dd0c0dSLois Curfman McInnes 954c8dd0c0dSLois Curfman McInnes Input parameters for func: 955c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 956c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 9578f99978dSLois Curfman McInnes - x - current candidate iterate 958c8dd0c0dSLois Curfman McInnes 959c8dd0c0dSLois Curfman McInnes Output parameters for func: 960c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 961c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 962c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 963c8dd0c0dSLois Curfman McInnes 964c8dd0c0dSLois Curfman McInnes Level: advanced 965c8dd0c0dSLois Curfman McInnes 9668f99978dSLois Curfman McInnes Notes: 967b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 968b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 969b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 970b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 971ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 9728f99978dSLois Curfman McInnes 973b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 974b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 975b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 976ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 977ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 978ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 979ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 980b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 9818f99978dSLois Curfman McInnes 982c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 983c8dd0c0dSLois Curfman McInnes 984c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 985c8dd0c0dSLois Curfman McInnes @*/ 98663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 987c8dd0c0dSLois Curfman McInnes { 9886849ba73SBarry Smith PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*); 989c8dd0c0dSLois Curfman McInnes 990c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 991c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 992c8dd0c0dSLois Curfman McInnes if (f) { 993c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 994c8dd0c0dSLois Curfman McInnes } 995c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 996c8dd0c0dSLois Curfman McInnes } 997c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 9986849ba73SBarry Smith typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/ 999c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 10004a2ae208SSatish Balay #undef __FUNCT__ 10014a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS" 100263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx) 1003c8dd0c0dSLois Curfman McInnes { 1004c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 10054b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->CheckStep = func; 10064b27c08aSLois Curfman McInnes ((SNES_LS *)(snes->data))->checkP = checkctx; 1007c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 1008c8dd0c0dSLois Curfman McInnes } 1009c8dd0c0dSLois Curfman McInnes EXTERN_C_END 1010c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 101104d965bbSLois Curfman McInnes /* 10124b27c08aSLois Curfman McInnes SNESPrintHelp_LS - Prints all options for the SNES_LS method. 1013329e7e40SMatthew Knepley 1014329e7e40SMatthew Knepley Input Parameter: 1015329e7e40SMatthew Knepley . snes - the SNES context 1016329e7e40SMatthew Knepley 1017329e7e40SMatthew Knepley Application Interface Routine: SNESPrintHelp() 1018329e7e40SMatthew Knepley */ 1019329e7e40SMatthew Knepley #undef __FUNCT__ 10204b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS" 10216849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p) 1022329e7e40SMatthew Knepley { 10234b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1024329e7e40SMatthew Knepley 1025329e7e40SMatthew Knepley PetscFunctionBegin; 10264b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 10274b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 10284b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 10294b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 10304b27c08aSLois Curfman McInnes (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 1031329e7e40SMatthew Knepley PetscFunctionReturn(0); 1032329e7e40SMatthew Knepley } 1033329e7e40SMatthew Knepley 1034329e7e40SMatthew Knepley /* 10354b27c08aSLois Curfman McInnes SNESView_LS - Prints info from the SNESLS data structure. 103604d965bbSLois Curfman McInnes 103704d965bbSLois Curfman McInnes Input Parameters: 103804d965bbSLois Curfman McInnes . SNES - the SNES context 103904d965bbSLois Curfman McInnes . viewer - visualization context 104004d965bbSLois Curfman McInnes 104104d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 104204d965bbSLois Curfman McInnes */ 10434a2ae208SSatish Balay #undef __FUNCT__ 10444b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS" 10456849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1046a935fc98SLois Curfman McInnes { 10474b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 10482fc52814SBarry Smith const char *cstr; 1049dfbe8321SBarry Smith PetscErrorCode ierr; 105032077d6dSBarry Smith PetscTruth iascii; 1051a935fc98SLois Curfman McInnes 10523a40ed3dSBarry Smith PetscFunctionBegin; 105332077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 105432077d6dSBarry Smith if (iascii) { 105519bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 105619bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 105719bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 105819bcc07fSBarry Smith else cstr = "unknown"; 1059b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1060b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 10615cd90555SBarry Smith } else { 106279a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 106319bcc07fSBarry Smith } 10643a40ed3dSBarry Smith PetscFunctionReturn(0); 1065a935fc98SLois Curfman McInnes } 106604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 106704d965bbSLois Curfman McInnes /* 10684b27c08aSLois Curfman McInnes SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 106904d965bbSLois Curfman McInnes 107004d965bbSLois Curfman McInnes Input Parameter: 107104d965bbSLois Curfman McInnes . snes - the SNES context 107204d965bbSLois Curfman McInnes 107304d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 107404d965bbSLois Curfman McInnes */ 10754a2ae208SSatish Balay #undef __FUNCT__ 10764b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS" 10776849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 10785e42470aSBarry Smith { 10794b27c08aSLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 1080e5999256SBarry Smith const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1081dfbe8321SBarry Smith PetscErrorCode ierr; 1082a7cc72afSBarry Smith PetscInt indx; 1083f1af5d2fSBarry Smith PetscTruth flg; 10845e42470aSBarry Smith 10853a40ed3dSBarry Smith PetscFunctionBegin; 1086b0a32e0cSBarry Smith ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 10874b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 10884b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 10894b27c08aSLois Curfman McInnes ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1090186905e3SBarry Smith 109122e36657SBarry Smith ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 109225cdf11fSBarry Smith if (flg) { 109322e36657SBarry Smith switch (indx) { 1094b49fd9e1SBarry Smith case 0: 1095af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1096b49fd9e1SBarry Smith break; 1097b49fd9e1SBarry Smith case 1: 1098af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1099b49fd9e1SBarry Smith break; 1100b49fd9e1SBarry Smith case 2: 1101af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1102b49fd9e1SBarry Smith break; 1103b49fd9e1SBarry Smith case 3: 1104af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1105b49fd9e1SBarry Smith break; 11065e42470aSBarry Smith } 11075e42470aSBarry Smith } 1108b0a32e0cSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 11093a40ed3dSBarry Smith PetscFunctionReturn(0); 11105e42470aSBarry Smith } 111104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 1112ebe3b25bSBarry Smith /*MC 1113ebe3b25bSBarry Smith SNESLS - Newton based nonlinear solver that uses a line search 111404d965bbSLois Curfman McInnes 1115ebe3b25bSBarry Smith Options Database: 1116ebe3b25bSBarry Smith + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1117ebe3b25bSBarry Smith . -snes_ls_alpha <alpha> - Sets alpha 1118ebe3b25bSBarry Smith . -snes_ls_maxstep <max> - Sets maxstep 1119ebe3b25bSBarry Smith - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1120ebe3b25bSBarry Smith will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1121ebe3b25bSBarry Smith the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 112204d965bbSLois Curfman McInnes 1123ebe3b25bSBarry Smith Notes: This is the default nonlinear solver in SNES 112404d965bbSLois Curfman McInnes 1125ee3001cbSBarry Smith Level: beginner 1126ee3001cbSBarry Smith 1127ebe3b25bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(), 1128ebe3b25bSBarry Smith SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(), 1129ebe3b25bSBarry Smith SNESSetLineSearch(), SNESNoLineSearchNoNorms() 1130ebe3b25bSBarry Smith 1131ebe3b25bSBarry Smith M*/ 1132fb2e594dSBarry Smith EXTERN_C_BEGIN 11334a2ae208SSatish Balay #undef __FUNCT__ 11344b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS" 113563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 11365e42470aSBarry Smith { 1137dfbe8321SBarry Smith PetscErrorCode ierr; 11384b27c08aSLois Curfman McInnes SNES_LS *neP; 11395e42470aSBarry Smith 11403a40ed3dSBarry Smith PetscFunctionBegin; 11414b27c08aSLois Curfman McInnes snes->setup = SNESSetUp_LS; 11424b27c08aSLois Curfman McInnes snes->solve = SNESSolve_LS; 11434b27c08aSLois Curfman McInnes snes->destroy = SNESDestroy_LS; 11444b27c08aSLois Curfman McInnes snes->converged = SNESConverged_LS; 11454b27c08aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_LS; 11464b27c08aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_LS; 11474b27c08aSLois Curfman McInnes snes->view = SNESView_LS; 11485baf8537SBarry Smith snes->nwork = 0; 11495e42470aSBarry Smith 11504b27c08aSLois Curfman McInnes ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 115152e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr); 11525e42470aSBarry Smith snes->data = (void*)neP; 11535e42470aSBarry Smith neP->alpha = 1.e-4; 11545e42470aSBarry Smith neP->maxstep = 1.e8; 11555e42470aSBarry Smith neP->steptol = 1.e-12; 11565e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1157c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1158c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1159c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 116082bf6240SBarry Smith 11615d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 11625d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 116382bf6240SBarry Smith 11643a40ed3dSBarry Smith PetscFunctionReturn(0); 11655e42470aSBarry Smith } 1166fb2e594dSBarry Smith EXTERN_C_END 116782bf6240SBarry Smith 116882bf6240SBarry Smith 116982bf6240SBarry Smith 1170