1*5d1a10b1SBarry Smith /*$Id: ls.c,v 1.158 2000/07/24 03:50:34 bsmith Exp bsmith $*/ 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 */ 118a5d9ee4SBarry Smith #undef __FUNC__ 128a5d9ee4SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckLocalMin_Private"></a>*/"SNESLSCheckLocalMin_Private" 138a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 148a5d9ee4SBarry Smith { 158a5d9ee4SBarry Smith PetscReal a1; 168a5d9ee4SBarry Smith int 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); 26*5d1a10b1SBarry Smith PLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm); 278a5d9ee4SBarry Smith if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 2836669109SBarry Smith } else { 2936669109SBarry Smith Vec work; 3036669109SBarry Smith Scalar 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); 40*5d1a10b1SBarry Smith PLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1); 4136669109SBarry Smith if (a1 < 1.e-4) *ismin = PETSC_TRUE; 4236669109SBarry Smith } 438a5d9ee4SBarry Smith PetscFunctionReturn(0); 448a5d9ee4SBarry Smith } 458a5d9ee4SBarry Smith 4604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------- 4704d965bbSLois Curfman McInnes 4804d965bbSLois Curfman McInnes This file implements a truncated Newton method with a line search, 4904d965bbSLois Curfman McInnes for solving a system of nonlinear equations, using the SLES, Vec, 5004d965bbSLois Curfman McInnes and Mat interfaces for linear solvers, vectors, and matrices, 5104d965bbSLois Curfman McInnes respectively. 5204d965bbSLois Curfman McInnes 53fe6c6bacSLois Curfman McInnes The following basic routines are required for each nonlinear solver: 5404d965bbSLois Curfman McInnes SNESCreate_XXX() - Creates a nonlinear solver context 5504d965bbSLois Curfman McInnes SNESSetFromOptions_XXX() - Sets runtime options 56fe6c6bacSLois Curfman McInnes SNESSolve_XXX() - Solves the nonlinear system 5704d965bbSLois Curfman McInnes SNESDestroy_XXX() - Destroys the nonlinear solver context 5804d965bbSLois Curfman McInnes The suffix "_XXX" denotes a particular implementation, in this case 5904d965bbSLois Curfman McInnes we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 6004d965bbSLois Curfman McInnes systems of nonlinear equations (EQ) with a line search (LS) method. 6104d965bbSLois Curfman McInnes These routines are actually called via the common user interface 6204d965bbSLois Curfman McInnes routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 6304d965bbSLois Curfman McInnes SNESDestroy(), so the application code interface remains identical 6404d965bbSLois Curfman McInnes for all nonlinear solvers. 6504d965bbSLois Curfman McInnes 6604d965bbSLois Curfman McInnes Another key routine is: 6704d965bbSLois Curfman McInnes SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 6804d965bbSLois Curfman McInnes by setting data structures and options. The interface routine SNESSetUp() 6904d965bbSLois Curfman McInnes is not usually called directly by the user, but instead is called by 7004d965bbSLois Curfman McInnes SNESSolve() if necessary. 7104d965bbSLois Curfman McInnes 7204d965bbSLois Curfman McInnes Additional basic routines are: 7304d965bbSLois Curfman McInnes SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 7404d965bbSLois Curfman McInnes SNESView_XXX() - Prints details of runtime options that 7504d965bbSLois Curfman McInnes have actually been used. 7604d965bbSLois Curfman McInnes These are called by application codes via the interface routines 7704d965bbSLois Curfman McInnes SNESPrintHelp() and SNESView(). 7804d965bbSLois Curfman McInnes 7904d965bbSLois Curfman McInnes The various types of solvers (preconditioners, Krylov subspace methods, 8004d965bbSLois Curfman McInnes nonlinear solvers, timesteppers) are all organized similarly, so the 8104d965bbSLois Curfman McInnes above description applies to these categories also. 8204d965bbSLois Curfman McInnes 8304d965bbSLois Curfman McInnes -------------------------------------------------------------------- */ 845e42470aSBarry Smith /* 8504d965bbSLois Curfman McInnes SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 8604d965bbSLois Curfman McInnes method with a line search. 875e76c431SBarry Smith 8804d965bbSLois Curfman McInnes Input Parameters: 8904d965bbSLois Curfman McInnes . snes - the SNES context 905e76c431SBarry Smith 9104d965bbSLois Curfman McInnes Output Parameter: 92c2a78f1aSLois Curfman McInnes . outits - number of iterations until termination 9304d965bbSLois Curfman McInnes 9404d965bbSLois Curfman McInnes Application Interface Routine: SNESSolve() 955e76c431SBarry Smith 965e76c431SBarry Smith Notes: 975e76c431SBarry Smith This implements essentially a truncated Newton method with a 985e76c431SBarry Smith line search. By default a cubic backtracking line search 995e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 1005e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 101393d2d9aSLois Curfman McInnes and Schnabel. 1025e76c431SBarry Smith */ 1035615d1e5SSatish Balay #undef __FUNC__ 104b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS" 105f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits) 1065e76c431SBarry Smith { 1076831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 10884c569c9SBarry Smith int maxits,i,ierr,lits,lsfail; 109112a2221SBarry Smith MatStructure flg = DIFFERENT_NONZERO_PATTERN; 110329f5518SBarry Smith PetscReal fnorm,gnorm,xnorm,ynorm; 1115e42470aSBarry Smith Vec Y,X,F,G,W,TMP; 1125e76c431SBarry Smith 1133a40ed3dSBarry Smith PetscFunctionBegin; 114184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 115184914b5SBarry Smith 1165e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 1175e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 11839e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 1195e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 1205e42470aSBarry Smith G = snes->work[1]; 1215e42470aSBarry Smith W = snes->work[2]; 1225e76c431SBarry Smith 1234c49b128SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1244c49b128SBarry Smith snes->iter = 0; 1254c49b128SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1265334005bSBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 127cddf8d76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1280f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1295e42470aSBarry Smith snes->norm = fnorm; 1300f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 13184c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,0); 13294a424c1SBarry Smith SNESMonitor(snes,0,fnorm); 1335e76c431SBarry Smith 134184914b5SBarry Smith if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 13594a9d846SBarry Smith 136d034289bSBarry Smith /* set parameter for default relative tolerance convergence test */ 137d034289bSBarry Smith snes->ttol = fnorm*snes->rtol; 138d034289bSBarry Smith 1395e76c431SBarry Smith for (i=0; i<maxits; i++) { 1405e76c431SBarry Smith 141ea4d3ed3SLois Curfman McInnes /* Solve J Y = F, where J is Jacobian matrix */ 1425334005bSBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 1435334005bSBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 14478b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 14590cbd584SBarry Smith /* should check what happened to the linear solve? */ 1463505fcc1SBarry Smith snes->linear_its += lits; 147af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 148ea4d3ed3SLois Curfman McInnes 149ea4d3ed3SLois Curfman McInnes /* Compute a (scaled) negative update in the line search routine: 150ea4d3ed3SLois Curfman McInnes Y <- X - lambda*Y 151ea4d3ed3SLois Curfman McInnes and evaluate G(Y) = function(Y)) 152ea4d3ed3SLois Curfman McInnes */ 15381b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 154af2d14d2SLois Curfman McInnes ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 155af1ccdebSLois Curfman McInnes PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 1563505fcc1SBarry Smith if (lsfail) { 1578a5d9ee4SBarry Smith PetscTruth ismin; 1583505fcc1SBarry Smith snes->nfailures++; 1593505fcc1SBarry Smith snes->reason = SNES_DIVERGED_LS_FAILURE; 1608a5d9ee4SBarry Smith ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 1618a5d9ee4SBarry Smith if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 1623505fcc1SBarry Smith break; 1633505fcc1SBarry Smith } 1645e76c431SBarry Smith 16539e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 16639e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 1675e76c431SBarry Smith fnorm = gnorm; 1685e76c431SBarry Smith 1690f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 17084c569c9SBarry Smith snes->iter = i+1; 1715e42470aSBarry Smith snes->norm = fnorm; 1720f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 17384c569c9SBarry Smith SNESLogConvHistory(snes,fnorm,lits); 17494a424c1SBarry Smith SNESMonitor(snes,i+1,fnorm); 1755e76c431SBarry Smith 1765e76c431SBarry Smith /* Test for convergence */ 17729e0b56fSBarry Smith if (snes->converged) { 17829e0b56fSBarry Smith ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 1793505fcc1SBarry Smith ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1803505fcc1SBarry Smith if (snes->reason) { 18116c95cacSBarry Smith break; 18216c95cacSBarry Smith } 18316c95cacSBarry Smith } 18429e0b56fSBarry Smith } 18539e2f89bSBarry Smith if (X != snes->vec_sol) { 186393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 187e15f7bb5SBarry Smith } 188e15f7bb5SBarry Smith if (F != snes->vec_func) { 189e15f7bb5SBarry Smith ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 190e15f7bb5SBarry Smith } 19139e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 19239e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 19352392280SLois Curfman McInnes if (i == maxits) { 194981c4779SBarry Smith PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 19552392280SLois Curfman McInnes i--; 1963505fcc1SBarry Smith snes->reason = SNES_DIVERGED_MAX_IT; 19752392280SLois Curfman McInnes } 1980f5bd95cSBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1990f5bd95cSBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2005e42470aSBarry Smith *outits = i+1; 2013a40ed3dSBarry Smith PetscFunctionReturn(0); 2025e76c431SBarry Smith } 20304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 20404d965bbSLois Curfman McInnes /* 20504d965bbSLois Curfman McInnes SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 2066831982aSBarry Smith of the SNESEQLS nonlinear solver. 20704d965bbSLois Curfman McInnes 20804d965bbSLois Curfman McInnes Input Parameter: 20904d965bbSLois Curfman McInnes . snes - the SNES context 21004d965bbSLois Curfman McInnes . x - the solution vector 21104d965bbSLois Curfman McInnes 21204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetUp() 21304d965bbSLois Curfman McInnes 21404d965bbSLois Curfman McInnes Notes: 21504d965bbSLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 21604d965bbSLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 21704d965bbSLois Curfman McInnes the call to SNESSolve(). 21804d965bbSLois Curfman McInnes */ 2195615d1e5SSatish Balay #undef __FUNC__ 220b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 221f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes) 2225e76c431SBarry Smith { 2235e42470aSBarry Smith int ierr; 2243a40ed3dSBarry Smith 2253a40ed3dSBarry Smith PetscFunctionBegin; 22681b6cf68SLois Curfman McInnes snes->nwork = 4; 227d7e8b826SBarry Smith ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 2285e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 22981b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 2303a40ed3dSBarry Smith PetscFunctionReturn(0); 2315e76c431SBarry Smith } 23204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 23304d965bbSLois Curfman McInnes /* 2346831982aSBarry Smith SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 23504d965bbSLois Curfman McInnes with SNESCreate_EQ_LS(). 23604d965bbSLois Curfman McInnes 23704d965bbSLois Curfman McInnes Input Parameter: 23804d965bbSLois Curfman McInnes . snes - the SNES context 23904d965bbSLois Curfman McInnes 240de49b36dSLois Curfman McInnes Application Interface Routine: SNESDestroy() 24104d965bbSLois Curfman McInnes */ 2425615d1e5SSatish Balay #undef __FUNC__ 243b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 244e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes) 2455e76c431SBarry Smith { 246393d2d9aSLois Curfman McInnes int ierr; 2473a40ed3dSBarry Smith 2483a40ed3dSBarry Smith PetscFunctionBegin; 2495baf8537SBarry Smith if (snes->nwork) { 2504b0e389bSBarry Smith ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 25121c89e3eSBarry Smith } 2525c704376SSatish Balay ierr = PetscFree(snes->data);CHKERRQ(ierr); 2533a40ed3dSBarry Smith PetscFunctionReturn(0); 2545e76c431SBarry Smith } 25504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 2565615d1e5SSatish Balay #undef __FUNC__ 257b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 25882bf6240SBarry Smith 2594b828684SBarry Smith /*@C 2605e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 2615e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 2625e76c431SBarry Smith to serve as a template and is not recommended for general use. 2635e76c431SBarry Smith 264c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 265c7afd0dbSLois Curfman McInnes 2665e76c431SBarry Smith Input Parameters: 267c7afd0dbSLois Curfman McInnes + snes - nonlinear context 268af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 2695e76c431SBarry Smith . x - current iterate 2705e76c431SBarry Smith . f - residual evaluated at x 2715e76c431SBarry Smith . y - search direction (contains new iterate on output) 2725e76c431SBarry Smith . w - work vector 273c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 2745e76c431SBarry Smith 275c4a48953SLois Curfman McInnes Output Parameters: 276c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 2775e76c431SBarry Smith . y - new iterate (contains search direction on input) 2785e76c431SBarry Smith . gnorm - 2-norm of g 2795e76c431SBarry Smith . ynorm - 2-norm of search length 280c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 281fee21e36SBarry Smith 282c4a48953SLois Curfman McInnes Options Database Key: 283c7afd0dbSLois Curfman McInnes . -snes_eq_ls basic - Activates SNESNoLineSearch() 284c4a48953SLois Curfman McInnes 28536851e7fSLois Curfman McInnes Level: advanced 28636851e7fSLois Curfman McInnes 28728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 28828ae5a21SLois Curfman McInnes 289f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 290af2d14d2SLois Curfman McInnes SNESSetLineSearch(), SNESNoLineSearchNoNorms() 2915e76c431SBarry Smith @*/ 292*5d1a10b1SBarry Smith int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 2935e76c431SBarry Smith { 2945e42470aSBarry Smith int ierr; 2955334005bSBarry Smith Scalar mone = -1.0; 2966831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 2978f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 2985334005bSBarry Smith 2993a40ed3dSBarry Smith PetscFunctionBegin; 300761aaf1bSLois Curfman McInnes *flag = 0; 3017857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 302a3b38805SLois Curfman McInnes ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 303ea4d3ed3SLois Curfman McInnes ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3048f99978dSLois Curfman McInnes if (neP->CheckStep) { 3058f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3068f99978dSLois Curfman McInnes } 307ea4d3ed3SLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 308a3b38805SLois Curfman McInnes ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 3097857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3103a40ed3dSBarry Smith PetscFunctionReturn(0); 3115e76c431SBarry Smith } 31204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 3135615d1e5SSatish Balay #undef __FUNC__ 314b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 31582bf6240SBarry Smith 31629e0b56fSBarry Smith /*@C 31729e0b56fSBarry Smith SNESNoLineSearchNoNorms - This routine is not a line search at 31829e0b56fSBarry Smith all; it simply uses the full Newton step. This version does not 31929e0b56fSBarry Smith even compute the norm of the function or search direction; this 32029e0b56fSBarry Smith is intended only when you know the full step is fine and are 32129e0b56fSBarry Smith not checking for convergence of the nonlinear iteration (for 322c7afd0dbSLois Curfman McInnes example, you are running always for a fixed number of Newton steps). 323c7afd0dbSLois Curfman McInnes 324c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 32529e0b56fSBarry Smith 32629e0b56fSBarry Smith Input Parameters: 327c7afd0dbSLois Curfman McInnes + snes - nonlinear context 328af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 32929e0b56fSBarry Smith . x - current iterate 33029e0b56fSBarry Smith . f - residual evaluated at x 33129e0b56fSBarry Smith . y - search direction (contains new iterate on output) 33229e0b56fSBarry Smith . w - work vector 333c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 33429e0b56fSBarry Smith 33529e0b56fSBarry Smith Output Parameters: 336c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 33729e0b56fSBarry Smith . gnorm - not changed 33829e0b56fSBarry Smith . ynorm - not changed 339c7afd0dbSLois Curfman McInnes - flag - set to 0, indicating a successful line search 340fee21e36SBarry Smith 34129e0b56fSBarry Smith Options Database Key: 342c7afd0dbSLois Curfman McInnes . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 34329e0b56fSBarry Smith 3448cbba510SBarry Smith Notes: 345ea56f5baSLois Curfman McInnes SNESNoLineSearchNoNorms() must be used in conjunction with 346ea56f5baSLois Curfman McInnes either the options 347ea56f5baSLois Curfman McInnes $ -snes_no_convergence_test -snes_max_it <its> 348ea56f5baSLois Curfman McInnes or alternatively a user-defined custom test set via 349329f5518SBarry Smith SNESSetConvergenceTest(); or a -snes_max_it of 1, 350329f5518SBarry Smith otherwise, the SNES solver will generate an error. 351329f5518SBarry Smith 352329f5518SBarry Smith During the final iteration this will not evaluate the function at 353329f5518SBarry Smith the solution point. This is to save a function evaluation while 354329f5518SBarry Smith using pseudo-timestepping. 3558cbba510SBarry Smith 356ea56f5baSLois Curfman McInnes The residual norms printed by monitoring routines such as 357ea56f5baSLois Curfman McInnes SNESDefaultMonitor() (as activated via -snes_monitor) will not be 358ea56f5baSLois Curfman McInnes correct, since they are not computed. 359ea56f5baSLois Curfman McInnes 360ea56f5baSLois Curfman McInnes Level: advanced 3618cbba510SBarry Smith 36229e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic 36329e0b56fSBarry Smith 36429e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 36529e0b56fSBarry Smith SNESSetLineSearch(), SNESNoLineSearch() 36629e0b56fSBarry Smith @*/ 367*5d1a10b1SBarry Smith int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 36829e0b56fSBarry Smith { 36929e0b56fSBarry Smith int ierr; 37029e0b56fSBarry Smith Scalar mone = -1.0; 3716831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 3728f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 37329e0b56fSBarry Smith 3743a40ed3dSBarry Smith PetscFunctionBegin; 3758cbba510SBarry Smith *flag = 0; 37629e0b56fSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 37729e0b56fSBarry Smith ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 3788f99978dSLois Curfman McInnes if (neP->CheckStep) { 3798f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 3808f99978dSLois Curfman McInnes } 381329f5518SBarry Smith 382329f5518SBarry Smith /* don't evaluate function the last time through */ 383329f5518SBarry Smith if (snes->iter < snes->max_its-1) { 38429e0b56fSBarry Smith ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 385329f5518SBarry Smith } 38629e0b56fSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3873a40ed3dSBarry Smith PetscFunctionReturn(0); 38829e0b56fSBarry Smith } 38904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 39029e0b56fSBarry Smith #undef __FUNC__ 391b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 3924b828684SBarry Smith /*@C 393f525115eSLois Curfman McInnes SNESCubicLineSearch - Performs a cubic line search (default line search method). 3945e76c431SBarry Smith 395c7afd0dbSLois Curfman McInnes Collective on SNES 396c7afd0dbSLois Curfman McInnes 3975e76c431SBarry Smith Input Parameters: 398c7afd0dbSLois Curfman McInnes + snes - nonlinear context 399af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 4005e76c431SBarry Smith . x - current iterate 4015e76c431SBarry Smith . f - residual evaluated at x 4025e76c431SBarry Smith . y - search direction (contains new iterate on output) 4035e76c431SBarry Smith . w - work vector 404c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 4055e76c431SBarry Smith 406393d2d9aSLois Curfman McInnes Output Parameters: 407c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 4085e76c431SBarry Smith . y - new iterate (contains search direction on input) 4095e76c431SBarry Smith . gnorm - 2-norm of g 4105e76c431SBarry Smith . ynorm - 2-norm of search length 411c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 412fee21e36SBarry Smith 413c4a48953SLois Curfman McInnes Options Database Key: 414c7afd0dbSLois Curfman McInnes $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 415c4a48953SLois Curfman McInnes 4165e76c431SBarry Smith Notes: 4175e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 4185e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 4195e76c431SBarry Smith 42036851e7fSLois Curfman McInnes Level: advanced 42136851e7fSLois Curfman McInnes 42228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 42328ae5a21SLois Curfman McInnes 424af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 4255e76c431SBarry Smith @*/ 426*5d1a10b1SBarry Smith int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 4275e76c431SBarry Smith { 428406556e6SLois Curfman McInnes /* 429406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 430406556e6SLois Curfman McInnes minimization problem: 431406556e6SLois Curfman McInnes min z(x): R^n -> R, 432406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 433406556e6SLois Curfman McInnes */ 434406556e6SLois Curfman McInnes 435329f5518SBarry Smith PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 436329f5518SBarry Smith PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 437aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4385e42470aSBarry Smith Scalar cinitslope,clambda; 4396b5873e3SBarry Smith #endif 4405e42470aSBarry Smith int ierr,count; 4416831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 4425334005bSBarry Smith Scalar mone = -1.0,scale; 4438f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 4445e76c431SBarry Smith 4453a40ed3dSBarry Smith PetscFunctionBegin; 4467857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 447761aaf1bSLois Curfman McInnes *flag = 0; 4485e76c431SBarry Smith alpha = neP->alpha; 4495e76c431SBarry Smith maxstep = neP->maxstep; 4505e76c431SBarry Smith steptol = neP->steptol; 4515e76c431SBarry Smith 452cddf8d76SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 453a1c2b6eeSLois Curfman McInnes if (*ynorm < snes->atol) { 454a1c2b6eeSLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 455a1c2b6eeSLois Curfman McInnes *gnorm = fnorm; 456a1c2b6eeSLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 457a1c2b6eeSLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 458ad922ac9SBarry Smith goto theend1; 45994a9d846SBarry Smith } 4605e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 4615e42470aSBarry Smith scale = maxstep/(*ynorm); 462aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 46390cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 4646b5873e3SBarry Smith #else 46590cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 4666b5873e3SBarry Smith #endif 467393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 4685e76c431SBarry Smith *ynorm = maxstep; 4695e76c431SBarry Smith } 4705e76c431SBarry Smith minlambda = steptol/(*ynorm); 471a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 472aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 473a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 474329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 4755e42470aSBarry Smith #else 476a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 4775e42470aSBarry Smith #endif 4785e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 4795e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 4805e76c431SBarry Smith 481393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 4825334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 48378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 484313b5538SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 485*5d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 486393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 48794a424c1SBarry Smith PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 488ad922ac9SBarry Smith goto theend1; 4895e76c431SBarry Smith } 4905e76c431SBarry Smith 4915e76c431SBarry Smith /* Fit points with quadratic */ 492313b5538SBarry Smith lambda = 1.0; 493a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 4945e76c431SBarry Smith lambdaprev = lambda; 4955e76c431SBarry Smith gnormprev = *gnorm; 496329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 497ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 498ddd12767SBarry Smith else lambda = lambdatemp; 499393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 500ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 501aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 502ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5035e42470aSBarry Smith #else 504ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5055e42470aSBarry Smith #endif 50678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 507cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 508*5d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 509393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 510ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 511ad922ac9SBarry Smith goto theend1; 5125e76c431SBarry Smith } 5135e76c431SBarry Smith 5145e76c431SBarry Smith /* Fit points with cubic */ 5155e76c431SBarry Smith count = 1; 5165e76c431SBarry Smith while (1) { 5175e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 5182b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 51990cbd584SBarry Smith PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 520393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 521761aaf1bSLois Curfman McInnes *flag = -1; break; 5225e76c431SBarry Smith } 523*5d1a10b1SBarry Smith t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 524*5d1a10b1SBarry Smith t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 525ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5262b022350SLois Curfman McInnes b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 5275e76c431SBarry Smith d = b*b - 3*a*initslope; 5285e76c431SBarry Smith if (d < 0.0) d = 0.0; 5295e76c431SBarry Smith if (a == 0.0) { 5305e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 5315e76c431SBarry Smith } else { 5325e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 5335e76c431SBarry Smith } 5345e76c431SBarry Smith lambdaprev = lambda; 5355e76c431SBarry Smith gnormprev = *gnorm; 536329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 537329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 5385e76c431SBarry Smith else lambda = lambdatemp; 539393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 540ea4d3ed3SLois Curfman McInnes lambdaneg = -lambda; 541aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 542ea4d3ed3SLois Curfman McInnes clambda = lambdaneg; 543393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 5445e42470aSBarry Smith #else 545ea4d3ed3SLois Curfman McInnes ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 5465e42470aSBarry Smith #endif 54778b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 548cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 549*5d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 550393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 551ea4d3ed3SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 5522715565aSLois Curfman McInnes goto theend1; 5532b022350SLois Curfman McInnes } else { 5542b022350SLois Curfman McInnes PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 5555e76c431SBarry Smith } 5565e76c431SBarry Smith count++; 5575e76c431SBarry Smith } 558ad922ac9SBarry Smith theend1: 5598f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 5608f99978dSLois Curfman McInnes if (neP->CheckStep) { 5618f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 5628f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 5638f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 5648f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 5658f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 5668f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 5678f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 5688f99978dSLois Curfman McInnes } 5698f99978dSLois Curfman McInnes } 5707857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 5713a40ed3dSBarry Smith PetscFunctionReturn(0); 5725e76c431SBarry Smith } 57304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 5745615d1e5SSatish Balay #undef __FUNC__ 575b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 5764b828684SBarry Smith /*@C 577f525115eSLois Curfman McInnes SNESQuadraticLineSearch - Performs a quadratic line search. 5785e76c431SBarry Smith 579c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 580c7afd0dbSLois Curfman McInnes 5815e42470aSBarry Smith Input Parameters: 582c7afd0dbSLois Curfman McInnes + snes - the SNES context 583af2d14d2SLois Curfman McInnes . lsctx - optional context for line search (not used here) 5845e42470aSBarry Smith . x - current iterate 5855e42470aSBarry Smith . f - residual evaluated at x 5865e42470aSBarry Smith . y - search direction (contains new iterate on output) 5875e42470aSBarry Smith . w - work vector 588c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 5895e42470aSBarry Smith 590c4a48953SLois Curfman McInnes Output Parameters: 591c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 5925e42470aSBarry Smith . y - new iterate (contains search direction on input) 5935e42470aSBarry Smith . gnorm - 2-norm of g 5945e42470aSBarry Smith . ynorm - 2-norm of search length 595c7afd0dbSLois Curfman McInnes - flag - 0 if line search succeeds; -1 on failure. 596fee21e36SBarry Smith 597c4a48953SLois Curfman McInnes Options Database Key: 598c7afd0dbSLois Curfman McInnes . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 599c4a48953SLois Curfman McInnes 6005e42470aSBarry Smith Notes: 6016831982aSBarry Smith Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 6025e42470aSBarry Smith 60336851e7fSLois Curfman McInnes Level: advanced 60436851e7fSLois Curfman McInnes 605f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 606f59ffdedSLois Curfman McInnes 607af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 6085e42470aSBarry Smith @*/ 609*5d1a10b1SBarry Smith int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 6105e76c431SBarry Smith { 611406556e6SLois Curfman McInnes /* 612406556e6SLois Curfman McInnes Note that for line search purposes we work with with the related 613406556e6SLois Curfman McInnes minimization problem: 614406556e6SLois Curfman McInnes min z(x): R^n -> R, 615406556e6SLois Curfman McInnes where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 616406556e6SLois Curfman McInnes */ 617329f5518SBarry Smith PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 618aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6195e42470aSBarry Smith Scalar cinitslope,clambda; 6206b5873e3SBarry Smith #endif 6215e42470aSBarry Smith int ierr,count; 6226831982aSBarry Smith SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 6235334005bSBarry Smith Scalar mone = -1.0,scale; 6248f99978dSLois Curfman McInnes PetscTruth change_y = PETSC_FALSE; 6255e76c431SBarry Smith 6263a40ed3dSBarry Smith PetscFunctionBegin; 6277857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 628761aaf1bSLois Curfman McInnes *flag = 0; 6295e42470aSBarry Smith alpha = neP->alpha; 6305e42470aSBarry Smith maxstep = neP->maxstep; 6315e42470aSBarry Smith steptol = neP->steptol; 6325e76c431SBarry Smith 6333505fcc1SBarry Smith ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 634b37302c6SLois Curfman McInnes if (*ynorm < snes->atol) { 63594a9d846SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 636b37302c6SLois Curfman McInnes *gnorm = fnorm; 637b37302c6SLois Curfman McInnes ierr = VecCopy(x,y);CHKERRQ(ierr); 638b37302c6SLois Curfman McInnes ierr = VecCopy(f,g);CHKERRQ(ierr); 639ad922ac9SBarry Smith goto theend2; 64094a9d846SBarry Smith } 6415e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 6425e42470aSBarry Smith scale = maxstep/(*ynorm); 643393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y);CHKERRQ(ierr); 6445e42470aSBarry Smith *ynorm = maxstep; 6455e76c431SBarry Smith } 6465e42470aSBarry Smith minlambda = steptol/(*ynorm); 647a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 648aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 649a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 650329f5518SBarry Smith initslope = PetscRealPart(cinitslope); 6515e42470aSBarry Smith #else 652a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 6535e42470aSBarry Smith #endif 6545e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 6555e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 6565e42470aSBarry Smith 657393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w);CHKERRQ(ierr); 6585334005bSBarry Smith ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 65978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 660cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 661*5d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 662393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 66394a424c1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 664ad922ac9SBarry Smith goto theend2; 6655e42470aSBarry Smith } 6665e42470aSBarry Smith 6675e42470aSBarry Smith /* Fit points with quadratic */ 668313b5538SBarry Smith lambda = 1.0; 6695e42470aSBarry Smith count = 1; 6705e42470aSBarry Smith while (1) { 6715e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 672981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 673*5d1a10b1SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 674393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 675761aaf1bSLois Curfman McInnes *flag = -1; break; 6765e42470aSBarry Smith } 677a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 678329f5518SBarry Smith if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 679329f5518SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 680329f5518SBarry Smith else lambda = lambdatemp; 681393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w);CHKERRQ(ierr); 6823505fcc1SBarry Smith lambdaneg = -lambda; 683aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6843505fcc1SBarry Smith clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 6855e42470aSBarry Smith #else 6863505fcc1SBarry Smith ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 6875e42470aSBarry Smith #endif 68878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 689cddf8d76SBarry Smith ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 690*5d1a10b1SBarry Smith if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 691393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y);CHKERRQ(ierr); 692981c4779SBarry Smith PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 69306259719SBarry Smith break; 6945e42470aSBarry Smith } 6955e42470aSBarry Smith count++; 6965e42470aSBarry Smith } 697ad922ac9SBarry Smith theend2: 6988f99978dSLois Curfman McInnes /* Optional user-defined check for line search step validity */ 6998f99978dSLois Curfman McInnes if (neP->CheckStep) { 7008f99978dSLois Curfman McInnes ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 7018f99978dSLois Curfman McInnes if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 7028f99978dSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 7038f99978dSLois Curfman McInnes ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 7048f99978dSLois Curfman McInnes ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 7058f99978dSLois Curfman McInnes ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 7068f99978dSLois Curfman McInnes ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 7078f99978dSLois Curfman McInnes } 7088f99978dSLois Curfman McInnes } 7097857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 7103a40ed3dSBarry Smith PetscFunctionReturn(0); 7115e76c431SBarry Smith } 71204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 7135615d1e5SSatish Balay #undef __FUNC__ 714b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 715c9e6a524SLois Curfman McInnes /*@C 71677c4ece6SBarry Smith SNESSetLineSearch - Sets the line search routine to be used 7176831982aSBarry Smith by the method SNESEQLS. 7185e76c431SBarry Smith 7195e76c431SBarry Smith Input Parameters: 720c7afd0dbSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 721af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for use by line search 722c7afd0dbSLois Curfman McInnes - func - pointer to int function 7235e76c431SBarry Smith 724fee21e36SBarry Smith Collective on SNES 725fee21e36SBarry Smith 726c4a48953SLois Curfman McInnes Available Routines: 727c7afd0dbSLois Curfman McInnes + SNESCubicLineSearch() - default line search 728f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 729f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 730af2d14d2SLois Curfman McInnes - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 7315e76c431SBarry Smith 732c4a48953SLois Curfman McInnes Options Database Keys: 733af2d14d2SLois Curfman McInnes + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 734c7afd0dbSLois Curfman McInnes . -snes_eq_ls_alpha <alpha> - Sets alpha 735c7afd0dbSLois Curfman McInnes . -snes_eq_ls_maxstep <max> - Sets maxstep 736c7afd0dbSLois Curfman McInnes - -snes_eq_ls_steptol <steptol> - Sets steptol 737c4a48953SLois Curfman McInnes 7385e76c431SBarry Smith Calling sequence of func: 739bd208895SLois Curfman McInnes .vb 740af2d14d2SLois Curfman McInnes func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 741329f5518SBarry Smith PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 742bd208895SLois Curfman McInnes .ve 7435e76c431SBarry Smith 7445e76c431SBarry Smith Input parameters for func: 745c7afd0dbSLois Curfman McInnes + snes - nonlinear context 746af2d14d2SLois Curfman McInnes . lsctx - optional user-defined context for line search 7475e76c431SBarry Smith . x - current iterate 7485e76c431SBarry Smith . f - residual evaluated at x 7495e76c431SBarry Smith . y - search direction (contains new iterate on output) 7505e76c431SBarry Smith . w - work vector 751c7afd0dbSLois Curfman McInnes - fnorm - 2-norm of f 7525e76c431SBarry Smith 7535e76c431SBarry Smith Output parameters for func: 754c7afd0dbSLois Curfman McInnes + g - residual evaluated at new iterate y 7555e76c431SBarry Smith . y - new iterate (contains search direction on input) 7565e76c431SBarry Smith . gnorm - 2-norm of g 7575e76c431SBarry Smith . ynorm - 2-norm of search length 758c7afd0dbSLois Curfman McInnes - flag - set to 0 if the line search succeeds; a nonzero integer 759761aaf1bSLois Curfman McInnes on failure. 760f59ffdedSLois Curfman McInnes 76136851e7fSLois Curfman McInnes Level: advanced 76236851e7fSLois Curfman McInnes 763f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 764f59ffdedSLois Curfman McInnes 765*5d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 766*5d1a10b1SBarry Smith SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 7675e76c431SBarry Smith @*/ 768329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 7695e76c431SBarry Smith { 770329f5518SBarry Smith int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 77182bf6240SBarry Smith 7723a40ed3dSBarry Smith PetscFunctionBegin; 77394d884c6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 77482bf6240SBarry Smith if (f) { 775af2d14d2SLois Curfman McInnes ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 77682bf6240SBarry Smith } 7773a40ed3dSBarry Smith PetscFunctionReturn(0); 7785e76c431SBarry Smith } 77904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 780fb2e594dSBarry Smith EXTERN_C_BEGIN 78182bf6240SBarry Smith #undef __FUNC__ 782b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 783af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 784af2d14d2SLois Curfman McInnes double,double*,double*,int*),void *lsctx) 78582bf6240SBarry Smith { 78682bf6240SBarry Smith PetscFunctionBegin; 7876831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 7886831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 78982bf6240SBarry Smith PetscFunctionReturn(0); 79082bf6240SBarry Smith } 791fb2e594dSBarry Smith EXTERN_C_END 79204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 793c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 794b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 795c8dd0c0dSLois Curfman McInnes /*@C 796530e4296SLois Curfman McInnes SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 7976831982aSBarry Smith by the line search routine in the Newton-based method SNESEQLS. 798c8dd0c0dSLois Curfman McInnes 799c8dd0c0dSLois Curfman McInnes Input Parameters: 800c8dd0c0dSLois Curfman McInnes + snes - nonlinear context obtained from SNESCreate() 801c8dd0c0dSLois Curfman McInnes . func - pointer to int function 802c8dd0c0dSLois Curfman McInnes - checkctx - optional user-defined context for use by step checking routine 803c8dd0c0dSLois Curfman McInnes 804c8dd0c0dSLois Curfman McInnes Collective on SNES 805c8dd0c0dSLois Curfman McInnes 806c8dd0c0dSLois Curfman McInnes Calling sequence of func: 807c8dd0c0dSLois Curfman McInnes .vb 808b1ae27eaSLois Curfman McInnes int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 809c8dd0c0dSLois Curfman McInnes .ve 810b1ae27eaSLois Curfman McInnes where func returns an error code of 0 on success and a nonzero 811b1ae27eaSLois Curfman McInnes on failure. 812c8dd0c0dSLois Curfman McInnes 813c8dd0c0dSLois Curfman McInnes Input parameters for func: 814c8dd0c0dSLois Curfman McInnes + snes - nonlinear context 815c8dd0c0dSLois Curfman McInnes . checkctx - optional user-defined context for use by step checking routine 8168f99978dSLois Curfman McInnes - x - current candidate iterate 817c8dd0c0dSLois Curfman McInnes 818c8dd0c0dSLois Curfman McInnes Output parameters for func: 819c8dd0c0dSLois Curfman McInnes + x - current iterate (possibly modified) 820c8dd0c0dSLois Curfman McInnes - flag - flag indicating whether x has been modified (either 821c8dd0c0dSLois Curfman McInnes PETSC_TRUE of PETSC_FALSE) 822c8dd0c0dSLois Curfman McInnes 823c8dd0c0dSLois Curfman McInnes Level: advanced 824c8dd0c0dSLois Curfman McInnes 8258f99978dSLois Curfman McInnes Notes: 826b1ae27eaSLois Curfman McInnes SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 827b1ae27eaSLois Curfman McInnes iterate computed by the line search checking routine. In particular, 828b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 829b1ae27eaSLois Curfman McInnes to the checking routine, and then (3) compute the corresponding nonlinear 830ea56f5baSLois Curfman McInnes function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 8318f99978dSLois Curfman McInnes 832b1ae27eaSLois Curfman McInnes SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 833b1ae27eaSLois Curfman McInnes new iterate computed by the line search checking routine. In particular, 834b1ae27eaSLois Curfman McInnes these routines (1) compute a candidate iterate u_{i+1} as well as a 835ea56f5baSLois Curfman McInnes candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 836ea56f5baSLois Curfman McInnes routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 837ea56f5baSLois Curfman McInnes were made to the candidate iterate in the checking routine (as indicated 838ea56f5baSLois Curfman McInnes by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 839b1ae27eaSLois Curfman McInnes very costly, so use this feature with caution! 8408f99978dSLois Curfman McInnes 841c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine 842c8dd0c0dSLois Curfman McInnes 843c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch() 844c8dd0c0dSLois Curfman McInnes @*/ 8458f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 846c8dd0c0dSLois Curfman McInnes { 8478f99978dSLois Curfman McInnes int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 848c8dd0c0dSLois Curfman McInnes 849c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 850c8dd0c0dSLois Curfman McInnes ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 851c8dd0c0dSLois Curfman McInnes if (f) { 852c8dd0c0dSLois Curfman McInnes ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 853c8dd0c0dSLois Curfman McInnes } 854c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 855c8dd0c0dSLois Curfman McInnes } 856c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 857c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN 858c8dd0c0dSLois Curfman McInnes #undef __FUNC__ 859b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 8608f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 861c8dd0c0dSLois Curfman McInnes { 862c8dd0c0dSLois Curfman McInnes PetscFunctionBegin; 8636831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 8646831982aSBarry Smith ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 865c8dd0c0dSLois Curfman McInnes PetscFunctionReturn(0); 866c8dd0c0dSLois Curfman McInnes } 867c8dd0c0dSLois Curfman McInnes EXTERN_C_END 868c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 86904d965bbSLois Curfman McInnes /* 8706831982aSBarry Smith SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 87182bf6240SBarry Smith 87204d965bbSLois Curfman McInnes Input Parameter: 87304d965bbSLois Curfman McInnes . snes - the SNES context 87404d965bbSLois Curfman McInnes 87504d965bbSLois Curfman McInnes Application Interface Routine: SNESPrintHelp() 87604d965bbSLois Curfman McInnes */ 8775615d1e5SSatish Balay #undef __FUNC__ 878b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 879f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 8805e42470aSBarry Smith { 8816831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 882d132466eSBarry Smith int ierr; 8836daaf66cSBarry Smith 8843a40ed3dSBarry Smith PetscFunctionBegin; 8856831982aSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 886d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 887d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 888d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 889d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 8903a40ed3dSBarry Smith PetscFunctionReturn(0); 8915e42470aSBarry Smith } 89204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 89304d965bbSLois Curfman McInnes /* 8946831982aSBarry Smith SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 89504d965bbSLois Curfman McInnes 89604d965bbSLois Curfman McInnes Input Parameters: 89704d965bbSLois Curfman McInnes . SNES - the SNES context 89804d965bbSLois Curfman McInnes . viewer - visualization context 89904d965bbSLois Curfman McInnes 90004d965bbSLois Curfman McInnes Application Interface Routine: SNESView() 90104d965bbSLois Curfman McInnes */ 9025615d1e5SSatish Balay #undef __FUNC__ 903b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 904e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer) 905a935fc98SLois Curfman McInnes { 9066831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 90719bcc07fSBarry Smith char *cstr; 90851695f54SLois Curfman McInnes int ierr; 9096831982aSBarry Smith PetscTruth isascii; 910a935fc98SLois Curfman McInnes 9113a40ed3dSBarry Smith PetscFunctionBegin; 9126831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 9130f5bd95cSBarry Smith if (isascii) { 91419bcc07fSBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 91519bcc07fSBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 91619bcc07fSBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 91719bcc07fSBarry Smith else cstr = "unknown"; 918d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 919d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 9205cd90555SBarry Smith } else { 9210f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 92219bcc07fSBarry Smith } 9233a40ed3dSBarry Smith PetscFunctionReturn(0); 924a935fc98SLois Curfman McInnes } 92504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 92604d965bbSLois Curfman McInnes /* 9276831982aSBarry Smith SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 92804d965bbSLois Curfman McInnes 92904d965bbSLois Curfman McInnes Input Parameter: 93004d965bbSLois Curfman McInnes . snes - the SNES context 93104d965bbSLois Curfman McInnes 93204d965bbSLois Curfman McInnes Application Interface Routine: SNESSetFromOptions() 93304d965bbSLois Curfman McInnes */ 9345615d1e5SSatish Balay #undef __FUNC__ 935b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 936f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes) 9375e42470aSBarry Smith { 9386831982aSBarry Smith SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 9395e42470aSBarry Smith char ver[16]; 9405e42470aSBarry Smith double tmp; 941f1af5d2fSBarry Smith int ierr; 942f1af5d2fSBarry Smith PetscTruth flg; 9435e42470aSBarry Smith 9443a40ed3dSBarry Smith PetscFunctionBegin; 94509d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 94625cdf11fSBarry Smith if (flg) { 9475e42470aSBarry Smith ls->alpha = tmp; 9485e42470aSBarry Smith } 94909d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 95025cdf11fSBarry Smith if (flg) { 9515e42470aSBarry Smith ls->maxstep = tmp; 9525e42470aSBarry Smith } 95309d61ba7SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 95425cdf11fSBarry Smith if (flg) { 9555e42470aSBarry Smith ls->steptol = tmp; 9565e42470aSBarry Smith } 95709d61ba7SLois Curfman McInnes ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 95825cdf11fSBarry Smith if (flg) { 959c38d4ed2SBarry Smith PetscTruth isbasic,isnonorms,isquad,iscubic; 9600f5bd95cSBarry Smith 961c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 962c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 963c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 964c38d4ed2SBarry Smith ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 9650f5bd95cSBarry Smith 9660f5bd95cSBarry Smith if (isbasic) { 967af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 9680f5bd95cSBarry Smith } else if (isnonorms) { 969af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 9700f5bd95cSBarry Smith } else if (isquad) { 971af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 9720f5bd95cSBarry Smith } else if (iscubic) { 973af2d14d2SLois Curfman McInnes ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 9745e42470aSBarry Smith } 975a8c6a408SBarry Smith else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 9765e42470aSBarry Smith } 9773a40ed3dSBarry Smith PetscFunctionReturn(0); 9785e42470aSBarry Smith } 97904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */ 98004d965bbSLois Curfman McInnes /* 9816831982aSBarry Smith SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 9826831982aSBarry Smith SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 98304d965bbSLois Curfman McInnes context, SNES, that was created within SNESCreate(). 98404d965bbSLois Curfman McInnes 98504d965bbSLois Curfman McInnes 98604d965bbSLois Curfman McInnes Input Parameter: 98704d965bbSLois Curfman McInnes . snes - the SNES context 98804d965bbSLois Curfman McInnes 98904d965bbSLois Curfman McInnes Application Interface Routine: SNESCreate() 99004d965bbSLois Curfman McInnes */ 991fb2e594dSBarry Smith EXTERN_C_BEGIN 9925615d1e5SSatish Balay #undef __FUNC__ 993b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 994f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes) 9955e42470aSBarry Smith { 99682bf6240SBarry Smith int ierr; 9976831982aSBarry Smith SNES_EQ_LS *neP; 9985e42470aSBarry Smith 9993a40ed3dSBarry Smith PetscFunctionBegin; 1000a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1001a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1002a8c6a408SBarry Smith } 100382bf6240SBarry Smith 1004f63b844aSLois Curfman McInnes snes->setup = SNESSetUp_EQ_LS; 1005f63b844aSLois Curfman McInnes snes->solve = SNESSolve_EQ_LS; 1006f63b844aSLois Curfman McInnes snes->destroy = SNESDestroy_EQ_LS; 1007f63b844aSLois Curfman McInnes snes->converged = SNESConverged_EQ_LS; 1008f63b844aSLois Curfman McInnes snes->printhelp = SNESPrintHelp_EQ_LS; 1009f63b844aSLois Curfman McInnes snes->setfromoptions = SNESSetFromOptions_EQ_LS; 1010f63b844aSLois Curfman McInnes snes->view = SNESView_EQ_LS; 10115baf8537SBarry Smith snes->nwork = 0; 10125e42470aSBarry Smith 10136831982aSBarry Smith neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 10146831982aSBarry Smith PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 10155e42470aSBarry Smith snes->data = (void*)neP; 10165e42470aSBarry Smith neP->alpha = 1.e-4; 10175e42470aSBarry Smith neP->maxstep = 1.e8; 10185e42470aSBarry Smith neP->steptol = 1.e-12; 10195e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 1020c8dd0c0dSLois Curfman McInnes neP->lsP = PETSC_NULL; 1021c8dd0c0dSLois Curfman McInnes neP->CheckStep = PETSC_NULL; 1022c8dd0c0dSLois Curfman McInnes neP->checkP = PETSC_NULL; 102382bf6240SBarry Smith 1024*5d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 1025*5d1a10b1SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 102682bf6240SBarry Smith 10273a40ed3dSBarry Smith PetscFunctionReturn(0); 10285e42470aSBarry Smith } 1029fb2e594dSBarry Smith EXTERN_C_END 103082bf6240SBarry Smith 103182bf6240SBarry Smith 103282bf6240SBarry Smith 1033