15e76c431SBarry Smith #ifndef lint 2*a703fe33SLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.46 1995/10/12 04:20:05 bsmith Exp curfman $"; 35e76c431SBarry Smith #endif 45e76c431SBarry Smith 55e76c431SBarry Smith #include <math.h> 6fbe28522SBarry Smith #include "ls.h" 7b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 85e42470aSBarry Smith 95e42470aSBarry Smith /* 105e42470aSBarry Smith Implements a line search variant of Newton's Method 115e76c431SBarry Smith for solving systems of nonlinear equations. 125e76c431SBarry Smith 135e76c431SBarry Smith Input parameters: 145e42470aSBarry Smith . snes - nonlinear context obtained from SNESCreate() 155e76c431SBarry Smith 165e42470aSBarry Smith Output Parameters: 17a935fc98SLois Curfman McInnes . outits - Number of global iterations until termination. 185e76c431SBarry Smith 195e76c431SBarry Smith Notes: 205e76c431SBarry Smith This implements essentially a truncated Newton method with a 215e76c431SBarry Smith line search. By default a cubic backtracking line search 225e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 235e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 24393d2d9aSLois Curfman McInnes and Schnabel. 255e76c431SBarry Smith */ 265e76c431SBarry Smith 275e42470aSBarry Smith int SNESSolve_LS(SNES snes,int *outits) 285e76c431SBarry Smith { 295e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 30761aaf1bSLois Curfman McInnes int maxits, i, history_len, ierr, lits, lsfail; 31df60cc22SBarry Smith MatStructure flg = ALLMAT_DIFFERENT_NONZERO_PATTERN; 326b5873e3SBarry Smith double fnorm, gnorm, xnorm, ynorm, *history; 335e42470aSBarry Smith Vec Y, X, F, G, W, TMP; 345e76c431SBarry Smith 355e42470aSBarry Smith history = snes->conv_hist; /* convergence history */ 365e42470aSBarry Smith history_len = snes->conv_hist_len; /* convergence history length */ 375e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 385e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 3939e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 405e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 415e42470aSBarry Smith G = snes->work[1]; 425e42470aSBarry Smith W = snes->work[2]; 435e76c431SBarry Smith 4478b31e54SBarry Smith ierr = SNESComputeInitialGuess(snes,X); CHKERRQ(ierr); /* X <- X_0 */ 45393d2d9aSLois Curfman McInnes ierr = VecNorm(X,&xnorm); CHKERRQ(ierr); /* xnorm = || X || */ 4678b31e54SBarry Smith ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr); /* (+/-) F(X) */ 47393d2d9aSLois Curfman McInnes ierr = VecNorm(F,&fnorm); CHKERRQ(ierr); /* fnorm <- ||F|| */ 485e42470aSBarry Smith snes->norm = fnorm; 495e76c431SBarry Smith if (history && history_len > 0) history[0] = fnorm; 50ddd12767SBarry Smith if (snes->monitor) {ierr = (*snes->monitor)(snes,0,fnorm,snes->monP); CHKERRQ(ierr);} 515e76c431SBarry Smith 525e76c431SBarry Smith for ( i=0; i<maxits; i++ ) { 535e42470aSBarry Smith snes->iter = i+1; 545e76c431SBarry Smith 555e76c431SBarry Smith /* Solve J Y = -F, where J is Jacobian matrix */ 56ddd12767SBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); 57ddd12767SBarry Smith CHKERRQ(ierr); 58ddd12767SBarry Smith ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); 59ddd12767SBarry Smith CHKERRQ(ierr); 6078b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr); 6181b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr); 62ddd12767SBarry Smith ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 63761aaf1bSLois Curfman McInnes if (lsfail) snes->nfailures++; 645e76c431SBarry Smith 6539e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 6639e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 675e76c431SBarry Smith fnorm = gnorm; 685e76c431SBarry Smith 695e42470aSBarry Smith snes->norm = fnorm; 705e76c431SBarry Smith if (history && history_len > i+1) history[i+1] = fnorm; 71393d2d9aSLois Curfman McInnes ierr = VecNorm(X,&xnorm); CHKERRQ(ierr); /* xnorm = || X || */ 72ddd12767SBarry Smith if (snes->monitor) {ierr = (*snes->monitor)(snes,i+1,fnorm,snes->monP);CHKERRQ(ierr);} 735e76c431SBarry Smith 745e76c431SBarry Smith /* Test for convergence */ 75bbb6d6a8SBarry Smith if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) { 7639e2f89bSBarry Smith if (X != snes->vec_sol) { 77393d2d9aSLois Curfman McInnes ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr); 7839e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 7939e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 8039e2f89bSBarry Smith } 815e76c431SBarry Smith break; 825e76c431SBarry Smith } 835e76c431SBarry Smith } 8452392280SLois Curfman McInnes if (i == maxits) { 8552392280SLois Curfman McInnes PLogInfo((PetscObject)snes, 86ddd12767SBarry Smith "SNESSolve_LS: Maximum number of iterations has been reached: %d\n",maxits); 8752392280SLois Curfman McInnes i--; 8852392280SLois Curfman McInnes } 895e42470aSBarry Smith *outits = i+1; 905e42470aSBarry Smith return 0; 915e76c431SBarry Smith } 925e76c431SBarry Smith /* ------------------------------------------------------------ */ 935e42470aSBarry Smith int SNESSetUp_LS(SNES snes ) 945e76c431SBarry Smith { 955e42470aSBarry Smith int ierr; 9681b6cf68SLois Curfman McInnes snes->nwork = 4; 9778b31e54SBarry Smith ierr = VecGetVecs(snes->vec_sol,snes->nwork,&snes->work); CHKERRQ(ierr); 985e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work); 9981b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1005e42470aSBarry Smith return 0; 1015e76c431SBarry Smith } 1025e76c431SBarry Smith /* ------------------------------------------------------------ */ 1035e42470aSBarry Smith int SNESDestroy_LS(PetscObject obj) 1045e76c431SBarry Smith { 1055e42470aSBarry Smith SNES snes = (SNES) obj; 106393d2d9aSLois Curfman McInnes int ierr; 107393d2d9aSLois Curfman McInnes ierr = VecFreeVecs(snes->work,snes->nwork); CHKERRQ(ierr); 10878b31e54SBarry Smith PETSCFREE(snes->data); 1095e42470aSBarry Smith return 0; 1105e76c431SBarry Smith } 1115e76c431SBarry Smith /* ------------------------------------------------------------ */ 1125e76c431SBarry Smith /*ARGSUSED*/ 1134b828684SBarry Smith /*@C 1145e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 1155e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 1165e76c431SBarry Smith to serve as a template and is not recommended for general use. 1175e76c431SBarry Smith 1185e76c431SBarry Smith Input Parameters: 1195e42470aSBarry Smith . snes - nonlinear context 1205e76c431SBarry Smith . x - current iterate 1215e76c431SBarry Smith . f - residual evaluated at x 1225e76c431SBarry Smith . y - search direction (contains new iterate on output) 1235e76c431SBarry Smith . w - work vector 1245e76c431SBarry Smith . fnorm - 2-norm of f 1255e76c431SBarry Smith 126c4a48953SLois Curfman McInnes Output Parameters: 1275e76c431SBarry Smith . g - residual evaluated at new iterate y 1285e76c431SBarry Smith . y - new iterate (contains search direction on input) 1295e76c431SBarry Smith . gnorm - 2-norm of g 1305e76c431SBarry Smith . ynorm - 2-norm of search length 131761aaf1bSLois Curfman McInnes . flag - set to 0, indicating a successful line search 1325e76c431SBarry Smith 133c4a48953SLois Curfman McInnes Options Database Key: 134c4a48953SLois Curfman McInnes $ -snes_line_search basic 135c4a48953SLois Curfman McInnes 13628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 13728ae5a21SLois Curfman McInnes 138f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 139a935fc98SLois Curfman McInnes SNESSetLineSearchRoutine() 1405e76c431SBarry Smith @*/ 1415e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 142761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 1435e76c431SBarry Smith { 1445e42470aSBarry Smith int ierr; 1455e42470aSBarry Smith Scalar one = 1.0; 146761aaf1bSLois Curfman McInnes *flag = 0; 1477857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 148393d2d9aSLois Curfman McInnes ierr = VecNorm(y,ynorm); CHKERRQ(ierr); /* ynorm = || y || */ 149393d2d9aSLois Curfman McInnes ierr = VecAXPY(&one,x,y); CHKERRQ(ierr); /* y <- x + y */ 150393d2d9aSLois Curfman McInnes ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* (+/-) F(y) */ 151393d2d9aSLois Curfman McInnes ierr = VecNorm(g,gnorm); CHKERRQ(ierr); /* gnorm = || g || */ 1527857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 153761aaf1bSLois Curfman McInnes return 0; 1545e76c431SBarry Smith } 1555e76c431SBarry Smith /* ------------------------------------------------------------------ */ 1564b828684SBarry Smith /*@C 157f59ffdedSLois Curfman McInnes SNESCubicLineSearch - This routine performs a cubic line search and 158f59ffdedSLois Curfman McInnes is the default line search method. 1595e76c431SBarry Smith 1605e76c431SBarry Smith Input Parameters: 1615e42470aSBarry Smith . snes - nonlinear context 1625e76c431SBarry Smith . x - current iterate 1635e76c431SBarry Smith . f - residual evaluated at x 1645e76c431SBarry Smith . y - search direction (contains new iterate on output) 1655e76c431SBarry Smith . w - work vector 1665e76c431SBarry Smith . fnorm - 2-norm of f 1675e76c431SBarry Smith 168393d2d9aSLois Curfman McInnes Output Parameters: 1695e76c431SBarry Smith . g - residual evaluated at new iterate y 1705e76c431SBarry Smith . y - new iterate (contains search direction on input) 1715e76c431SBarry Smith . gnorm - 2-norm of g 1725e76c431SBarry Smith . ynorm - 2-norm of search length 173761aaf1bSLois Curfman McInnes . flag - 0 if line search succeeds; -1 on failure. 1745e76c431SBarry Smith 175c4a48953SLois Curfman McInnes Options Database Key: 176c4a48953SLois Curfman McInnes $ -snes_line_search cubic 177c4a48953SLois Curfman McInnes 1785e76c431SBarry Smith Notes: 1795e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 1805e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 1815e76c431SBarry Smith 18228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 18328ae5a21SLois Curfman McInnes 184f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine() 1855e76c431SBarry Smith @*/ 1865e42470aSBarry Smith int SNESCubicLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 187761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 1885e76c431SBarry Smith { 189ddd12767SBarry Smith double steptol, initslope,lambdaprev, gnormprev,a, b, d, t1, t2; 190ddd12767SBarry Smith double maxstep,minlambda,alpha,lambda,lambdatemp; 1916b5873e3SBarry Smith #if defined(PETSC_COMPLEX) 1925e42470aSBarry Smith Scalar cinitslope,clambda; 1936b5873e3SBarry Smith #endif 1945e42470aSBarry Smith int ierr, count; 1955e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 1965e42470aSBarry Smith Scalar one = 1.0,scale; 1975e76c431SBarry Smith 1987857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 199761aaf1bSLois Curfman McInnes *flag = 0; 2005e76c431SBarry Smith alpha = neP->alpha; 2015e76c431SBarry Smith maxstep = neP->maxstep; 2025e76c431SBarry Smith steptol = neP->steptol; 2035e76c431SBarry Smith 204393d2d9aSLois Curfman McInnes ierr = VecNorm(y,ynorm); CHKERRQ(ierr); 2055e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 2065e42470aSBarry Smith scale = maxstep/(*ynorm); 2076b5873e3SBarry Smith #if defined(PETSC_COMPLEX) 208ddd12767SBarry Smith PLogInfo((PetscObject)snes,"SNESCubicLineSearch: Scaling step by %g\n",real(scale)); 2096b5873e3SBarry Smith #else 210ddd12767SBarry Smith PLogInfo((PetscObject)snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 2116b5873e3SBarry Smith #endif 212393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 2135e76c431SBarry Smith *ynorm = maxstep; 2145e76c431SBarry Smith } 2155e76c431SBarry Smith minlambda = steptol/(*ynorm); 216*a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 2175e42470aSBarry Smith #if defined(PETSC_COMPLEX) 218*a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 2195e42470aSBarry Smith initslope = real(cinitslope); 2205e42470aSBarry Smith #else 221*a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 2225e42470aSBarry Smith #endif 2235e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 2245e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 2255e76c431SBarry Smith 226393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 227393d2d9aSLois Curfman McInnes ierr = VecAXPY(&one,x,w); CHKERRQ(ierr); 22878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 229393d2d9aSLois Curfman McInnes ierr = VecNorm(g,gnorm); 2305e76c431SBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* Sufficient reduction */ 231393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 2324b828684SBarry Smith PLogInfo((PetscObject)snes,"SNESCubicLineSearch: Using full step\n"); 233d93a2b8dSBarry Smith goto theend; 2345e76c431SBarry Smith } 2355e76c431SBarry Smith 2365e76c431SBarry Smith /* Fit points with quadratic */ 2375e76c431SBarry Smith lambda = 1.0; count = 0; 238*a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 2395e76c431SBarry Smith lambdaprev = lambda; 2405e76c431SBarry Smith gnormprev = *gnorm; 241ddd12767SBarry Smith if (lambdatemp <= .1*lambda) lambda = .1*lambda; 242ddd12767SBarry Smith else lambda = lambdatemp; 243393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 2445e42470aSBarry Smith #if defined(PETSC_COMPLEX) 245393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 2465e42470aSBarry Smith #else 247393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 2485e42470aSBarry Smith #endif 24978b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 250393d2d9aSLois Curfman McInnes ierr = VecNorm(g,gnorm); CHKERRQ(ierr); 2515e76c431SBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* sufficient reduction */ 252393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 2534b828684SBarry Smith PLogInfo((PetscObject)snes, 2544b828684SBarry Smith "SNESCubicLineSearch: Quadratically determined step, lambda %g\n",lambda); 255d93a2b8dSBarry Smith goto theend; 2565e76c431SBarry Smith } 2575e76c431SBarry Smith 2585e76c431SBarry Smith /* Fit points with cubic */ 2595e76c431SBarry Smith count = 1; 2605e76c431SBarry Smith while (1) { 2615e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 2624b828684SBarry Smith PLogInfo((PetscObject)snes, 2634b828684SBarry Smith "SNESCubicLineSearch:Unable to find good step length! %d \n",count); 2644b828684SBarry Smith PLogInfo((PetscObject)snes, 265ddd12767SBarry Smith "SNESCubicLineSearch:f %g fnew %g ynorm %g lambda %g\n",fnorm,*gnorm,*ynorm,lambda); 266393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 267761aaf1bSLois Curfman McInnes *flag = -1; break; 2685e76c431SBarry Smith } 2695e76c431SBarry Smith t1 = *gnorm - fnorm - lambda*initslope; 2705e76c431SBarry Smith t2 = gnormprev - fnorm - lambdaprev*initslope; 271ddd12767SBarry Smith a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 2725e76c431SBarry Smith b = (-lambdaprev*t1/(lambda*lambda) + 2735e76c431SBarry Smith lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 2745e76c431SBarry Smith d = b*b - 3*a*initslope; 2755e76c431SBarry Smith if (d < 0.0) d = 0.0; 2765e76c431SBarry Smith if (a == 0.0) { 2775e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 2785e76c431SBarry Smith } else { 2795e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 2805e76c431SBarry Smith } 2815e76c431SBarry Smith if (lambdatemp > .5*lambda) { 2825e76c431SBarry Smith lambdatemp = .5*lambda; 2835e76c431SBarry Smith } 2845e76c431SBarry Smith lambdaprev = lambda; 2855e76c431SBarry Smith gnormprev = *gnorm; 2865e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 2875e76c431SBarry Smith lambda = .1*lambda; 2885e76c431SBarry Smith } 2895e76c431SBarry Smith else lambda = lambdatemp; 290393d2d9aSLois Curfman McInnes ierr = VecCopy( x, w ); CHKERRQ(ierr); 2915e42470aSBarry Smith #if defined(PETSC_COMPLEX) 292393d2d9aSLois Curfman McInnes ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 2935e42470aSBarry Smith #else 294393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 2955e42470aSBarry Smith #endif 29678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 297393d2d9aSLois Curfman McInnes ierr = VecNorm(g,gnorm); CHKERRQ(ierr); 2985e76c431SBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* is reduction enough */ 299393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 3004b828684SBarry Smith PLogInfo((PetscObject)snes, 3014b828684SBarry Smith "SNESCubicLineSearch: Cubically determined step, lambda %g\n",lambda); 302761aaf1bSLois Curfman McInnes *flag = -1; break; 3035e76c431SBarry Smith } 3045e76c431SBarry Smith count++; 3055e76c431SBarry Smith } 306d93a2b8dSBarry Smith theend: 3077857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3085e42470aSBarry Smith return 0; 3095e76c431SBarry Smith } 31052392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 3114b828684SBarry Smith /*@C 3125e42470aSBarry Smith SNESQuadraticLineSearch - This routine performs a cubic line search. 3135e76c431SBarry Smith 3145e42470aSBarry Smith Input Parameters: 315f59ffdedSLois Curfman McInnes . snes - the SNES context 3165e42470aSBarry Smith . x - current iterate 3175e42470aSBarry Smith . f - residual evaluated at x 3185e42470aSBarry Smith . y - search direction (contains new iterate on output) 3195e42470aSBarry Smith . w - work vector 3205e42470aSBarry Smith . fnorm - 2-norm of f 3215e42470aSBarry Smith 322c4a48953SLois Curfman McInnes Output Parameters: 3235e42470aSBarry Smith . g - residual evaluated at new iterate y 3245e42470aSBarry Smith . y - new iterate (contains search direction on input) 3255e42470aSBarry Smith . gnorm - 2-norm of g 3265e42470aSBarry Smith . ynorm - 2-norm of search length 327761aaf1bSLois Curfman McInnes . flag - 0 if line search succeeds; -1 on failure. 3285e42470aSBarry Smith 329c4a48953SLois Curfman McInnes Options Database Key: 330c4a48953SLois Curfman McInnes $ -snes_line_search quadratic 331c4a48953SLois Curfman McInnes 3325e42470aSBarry Smith Notes: 333f59ffdedSLois Curfman McInnes Use SNESSetLineSearchRoutine() 334ddd12767SBarry Smith to set this routine within the SNES_EQ_NLS method. 3355e42470aSBarry Smith 336f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 337f59ffdedSLois Curfman McInnes 338f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine() 3395e42470aSBarry Smith @*/ 3405e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 341761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 3425e76c431SBarry Smith { 343ddd12767SBarry Smith double steptol,initslope,lambdaprev,gnormprev,maxstep,minlambda,alpha,lambda,lambdatemp; 3446b5873e3SBarry Smith #if defined(PETSC_COMPLEX) 3455e42470aSBarry Smith Scalar cinitslope,clambda; 3466b5873e3SBarry Smith #endif 3475e42470aSBarry Smith int ierr,count; 3485e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3495e42470aSBarry Smith Scalar one = 1.0,scale; 3505e76c431SBarry Smith 3517857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 352761aaf1bSLois Curfman McInnes *flag = 0; 3535e42470aSBarry Smith alpha = neP->alpha; 3545e42470aSBarry Smith maxstep = neP->maxstep; 3555e42470aSBarry Smith steptol = neP->steptol; 3565e76c431SBarry Smith 3575e42470aSBarry Smith VecNorm(y, ynorm ); 3585e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3595e42470aSBarry Smith scale = maxstep/(*ynorm); 360393d2d9aSLois Curfman McInnes ierr = VecScale(&scale,y); CHKERRQ(ierr); 3615e42470aSBarry Smith *ynorm = maxstep; 3625e76c431SBarry Smith } 3635e42470aSBarry Smith minlambda = steptol/(*ynorm); 364*a703fe33SLois Curfman McInnes ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr); 3655e42470aSBarry Smith #if defined(PETSC_COMPLEX) 366*a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr); 3675e42470aSBarry Smith initslope = real(cinitslope); 3685e42470aSBarry Smith #else 369*a703fe33SLois Curfman McInnes ierr = VecDot(f,w,&initslope); CHKERRQ(ierr); 3705e42470aSBarry Smith #endif 3715e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 3725e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 3735e42470aSBarry Smith 374393d2d9aSLois Curfman McInnes ierr = VecCopy(y,w); CHKERRQ(ierr); 375393d2d9aSLois Curfman McInnes ierr = VecAXPY(&one,x,w); CHKERRQ(ierr); 37678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 377393d2d9aSLois Curfman McInnes ierr = VecNorm(g,gnorm); CHKERRQ(ierr); 3785e42470aSBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* Sufficient reduction */ 379393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 3804b828684SBarry Smith PLogInfo((PetscObject)snes,"SNESQuadraticLineSearch: Using full step\n"); 381d93a2b8dSBarry Smith goto theend; 3825e42470aSBarry Smith } 3835e42470aSBarry Smith 3845e42470aSBarry Smith /* Fit points with quadratic */ 3855e42470aSBarry Smith lambda = 1.0; count = 0; 3865e42470aSBarry Smith count = 1; 3875e42470aSBarry Smith while (1) { 3885e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 3894b828684SBarry Smith PLogInfo((PetscObject)snes, 3904b828684SBarry Smith "SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 3914b828684SBarry Smith PLogInfo((PetscObject)snes, 392ddd12767SBarry Smith "SNESQuadraticLineSearch:f %g fnew %g ynorm %g lambda %g\n",fnorm,*gnorm,*ynorm,lambda); 393393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 394761aaf1bSLois Curfman McInnes *flag = -1; break; 3955e42470aSBarry Smith } 396*a703fe33SLois Curfman McInnes lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 3975e42470aSBarry Smith lambdaprev = lambda; 3985e42470aSBarry Smith gnormprev = *gnorm; 3995e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 4005e42470aSBarry Smith lambda = .1*lambda; 4015e42470aSBarry Smith } else lambda = lambdatemp; 402393d2d9aSLois Curfman McInnes ierr = VecCopy(x,w); CHKERRQ(ierr); 4035e42470aSBarry Smith #if defined(PETSC_COMPLEX) 404393d2d9aSLois Curfman McInnes clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr); 4055e42470aSBarry Smith #else 406393d2d9aSLois Curfman McInnes ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr); 4075e42470aSBarry Smith #endif 40878b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 409393d2d9aSLois Curfman McInnes ierr = VecNorm(g,gnorm); CHKERRQ(ierr); 4105e42470aSBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* sufficient reduction */ 411393d2d9aSLois Curfman McInnes ierr = VecCopy(w,y); CHKERRQ(ierr); 4124b828684SBarry Smith PLogInfo((PetscObject)snes, 413ddd12767SBarry Smith "SNESQuadraticLineSearch:Quadratically determined step, lambda %g\n",lambda); 41406259719SBarry Smith break; 4155e42470aSBarry Smith } 4165e42470aSBarry Smith count++; 4175e42470aSBarry Smith } 418d93a2b8dSBarry Smith theend: 4197857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4205e42470aSBarry Smith return 0; 4215e76c431SBarry Smith } 4225e76c431SBarry Smith /* ------------------------------------------------------------ */ 423b1f0a012SBarry Smith /*@ 4245e42470aSBarry Smith SNESSetLineSearchRoutine - Sets the line search routine to be used 425fbe28522SBarry Smith by the method SNES_LS. 4265e76c431SBarry Smith 4275e76c431SBarry Smith Input Parameters: 428eafb4bcbSBarry Smith . snes - nonlinear context obtained from SNESCreate() 4295e76c431SBarry Smith . func - pointer to int function 4305e76c431SBarry Smith 431c4a48953SLois Curfman McInnes Available Routines: 432f59ffdedSLois Curfman McInnes . SNESCubicLineSearch() - default line search 433f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 434f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 4355e76c431SBarry Smith 436c4a48953SLois Curfman McInnes Options Database Keys: 437c4a48953SLois Curfman McInnes $ -snes_line_search [basic,quadratic,cubic] 438c4a48953SLois Curfman McInnes 4395e76c431SBarry Smith Calling sequence of func: 440f59ffdedSLois Curfman McInnes func (SNES snes, Vec x, Vec f, Vec g, Vec y, 441761aaf1bSLois Curfman McInnes Vec w, double fnorm, double *ynorm, 442761aaf1bSLois Curfman McInnes double *gnorm, *flag) 4435e76c431SBarry Smith 4445e76c431SBarry Smith Input parameters for func: 4455e42470aSBarry Smith . snes - nonlinear context 4465e76c431SBarry Smith . x - current iterate 4475e76c431SBarry Smith . f - residual evaluated at x 4485e76c431SBarry Smith . y - search direction (contains new iterate on output) 4495e76c431SBarry Smith . w - work vector 4505e76c431SBarry Smith . fnorm - 2-norm of f 4515e76c431SBarry Smith 4525e76c431SBarry Smith Output parameters for func: 4535e76c431SBarry Smith . g - residual evaluated at new iterate y 4545e76c431SBarry Smith . y - new iterate (contains search direction on input) 4555e76c431SBarry Smith . gnorm - 2-norm of g 4565e76c431SBarry Smith . ynorm - 2-norm of search length 457761aaf1bSLois Curfman McInnes . flag - set to 0 if the line search succeeds; a nonzero integer 458761aaf1bSLois Curfman McInnes on failure. 459f59ffdedSLois Curfman McInnes 460f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 461f59ffdedSLois Curfman McInnes 462f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch() 4635e76c431SBarry Smith @*/ 4645e42470aSBarry Smith int SNESSetLineSearchRoutine(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec, 465761aaf1bSLois Curfman McInnes double,double *,double*,int*) ) 4665e76c431SBarry Smith { 467ddd12767SBarry Smith if ((snes)->type == SNES_EQ_NLS) 4685e42470aSBarry Smith ((SNES_LS *)(snes->data))->LineSearch = func; 4695e42470aSBarry Smith return 0; 4705e76c431SBarry Smith } 47152392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 4725e42470aSBarry Smith static int SNESPrintHelp_LS(SNES snes) 4735e42470aSBarry Smith { 4745e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 47552392280SLois Curfman McInnes char *p; 47652392280SLois Curfman McInnes if (snes->prefix) p = snes->prefix; else p = "-"; 47752392280SLois Curfman McInnes MPIU_printf(snes->comm," method ls (system of nonlinear equations):\n"); 47852392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search [basic,quadratic,cubic]\n",p); 47952392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search_alpha alpha (default %g)\n",p,ls->alpha); 48052392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search_maxstep max (default %g)\n",p,ls->maxstep); 48152392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search_steptol tol (default %g)\n",p,ls->steptol); 4826b5873e3SBarry Smith return 0; 4835e42470aSBarry Smith } 48452392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 485a935fc98SLois Curfman McInnes static int SNESView_LS(PetscObject obj,Viewer viewer) 486a935fc98SLois Curfman McInnes { 487a935fc98SLois Curfman McInnes SNES snes = (SNES)obj; 488a935fc98SLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 489d636dbe3SBarry Smith FILE *fd; 490a935fc98SLois Curfman McInnes char *cstring; 49151695f54SLois Curfman McInnes int ierr; 492a935fc98SLois Curfman McInnes 493a447f0b5SLois Curfman McInnes ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 494ddd12767SBarry Smith if (ls->LineSearch == SNESNoLineSearch) cstring = "SNESNoLineSearch"; 495ddd12767SBarry Smith else if (ls->LineSearch == SNESQuadraticLineSearch) cstring = "SNESQuadraticLineSearch"; 496ddd12767SBarry Smith else if (ls->LineSearch == SNESCubicLineSearch) cstring = "SNESCubicLineSearch"; 497ddd12767SBarry Smith else cstring = "unknown"; 498a935fc98SLois Curfman McInnes MPIU_fprintf(snes->comm,fd," line search variant: %s\n",cstring); 499a935fc98SLois Curfman McInnes MPIU_fprintf(snes->comm,fd," alpha=%g, maxstep=%g, steptol=%g\n", 500a935fc98SLois Curfman McInnes ls->alpha,ls->maxstep,ls->steptol); 501a935fc98SLois Curfman McInnes return 0; 502a935fc98SLois Curfman McInnes } 50352392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 5045e42470aSBarry Smith static int SNESSetFromOptions_LS(SNES snes) 5055e42470aSBarry Smith { 5065e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 5075e42470aSBarry Smith char ver[16]; 5085e42470aSBarry Smith double tmp; 5095e42470aSBarry Smith 510df60cc22SBarry Smith if (OptionsGetDouble(snes->prefix,"-snes_line_search_alpa",&tmp)) { 5115e42470aSBarry Smith ls->alpha = tmp; 5125e42470aSBarry Smith } 513df60cc22SBarry Smith if (OptionsGetDouble(snes->prefix,"-snes_line_search_maxstep",&tmp)) { 5145e42470aSBarry Smith ls->maxstep = tmp; 5155e42470aSBarry Smith } 516df60cc22SBarry Smith if (OptionsGetDouble(snes->prefix,"-snes_line_search_steptol",&tmp)) { 5175e42470aSBarry Smith ls->steptol = tmp; 5185e42470aSBarry Smith } 519df60cc22SBarry Smith if (OptionsGetString(snes->prefix,"-snes_line_search",ver,16)) { 52048d91487SBarry Smith if (!PetscStrcmp(ver,"basic")) { 5215e42470aSBarry Smith SNESSetLineSearchRoutine(snes,SNESNoLineSearch); 5225e42470aSBarry Smith } 52348d91487SBarry Smith else if (!PetscStrcmp(ver,"quadratic")) { 5245e42470aSBarry Smith SNESSetLineSearchRoutine(snes,SNESQuadraticLineSearch); 5255e42470aSBarry Smith } 52648d91487SBarry Smith else if (!PetscStrcmp(ver,"cubic")) { 5275e42470aSBarry Smith SNESSetLineSearchRoutine(snes,SNESCubicLineSearch); 5285e42470aSBarry Smith } 5298c48e012SBarry Smith else {SETERRQ(1,"SNESSetFromOptions_LS:Unknown line search");} 5305e42470aSBarry Smith } 5315e42470aSBarry Smith return 0; 5325e42470aSBarry Smith } 5335e42470aSBarry Smith /* ------------------------------------------------------------ */ 5345e42470aSBarry Smith int SNESCreate_LS(SNES snes ) 5355e42470aSBarry Smith { 5365e42470aSBarry Smith SNES_LS *neP; 5375e42470aSBarry Smith 538ddd12767SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 539ddd12767SBarry Smith SETERRQ(1,"SNESCreate_LS:For SNES_NONLINEAR_EQUATIONS only"); 54025c7317fSLois Curfman McInnes snes->type = SNES_EQ_NLS; 54149e3953dSBarry Smith snes->setup = SNESSetUp_LS; 54249e3953dSBarry Smith snes->solve = SNESSolve_LS; 5435e42470aSBarry Smith snes->destroy = SNESDestroy_LS; 544bbb6d6a8SBarry Smith snes->converged = SNESDefaultConverged; 54549e3953dSBarry Smith snes->printhelp = SNESPrintHelp_LS; 54649e3953dSBarry Smith snes->setfromoptions = SNESSetFromOptions_LS; 547a935fc98SLois Curfman McInnes snes->view = SNESView_LS; 5485e42470aSBarry Smith 54978b31e54SBarry Smith neP = PETSCNEW(SNES_LS); CHKPTRQ(neP); 550ff782a9fSLois Curfman McInnes PLogObjectMemory(snes,sizeof(SNES_LS)); 5515e42470aSBarry Smith snes->data = (void *) neP; 5525e42470aSBarry Smith neP->alpha = 1.e-4; 5535e42470aSBarry Smith neP->maxstep = 1.e8; 5545e42470aSBarry Smith neP->steptol = 1.e-12; 5555e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 5565e42470aSBarry Smith return 0; 5575e42470aSBarry Smith } 5585e42470aSBarry Smith 5595e42470aSBarry Smith 5605e42470aSBarry Smith 5615e42470aSBarry Smith 562