15e76c431SBarry Smith #ifndef lint 2*52392280SLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.34 1995/07/22 19:46:43 curfman Exp curfman $"; 35e76c431SBarry Smith #endif 45e76c431SBarry Smith 55e76c431SBarry Smith #include <math.h> 6fbe28522SBarry Smith #include "ls.h" 7a935fc98SLois Curfman McInnes #include "pviewer.h" 8bbb6d6a8SBarry Smith #if defined(HAVE_STRING_H) 9bbb6d6a8SBarry Smith #include <string.h> 10bbb6d6a8SBarry Smith #endif 115e42470aSBarry Smith 125e42470aSBarry Smith /* 135e42470aSBarry Smith Implements a line search variant of Newton's Method 145e76c431SBarry Smith for solving systems of nonlinear equations. 155e76c431SBarry Smith 165e76c431SBarry Smith Input parameters: 175e42470aSBarry Smith . snes - nonlinear context obtained from SNESCreate() 185e76c431SBarry Smith 195e42470aSBarry Smith Output Parameters: 20a935fc98SLois Curfman McInnes . outits - Number of global iterations until termination. 215e76c431SBarry Smith 225e76c431SBarry Smith Notes: 235e76c431SBarry Smith This implements essentially a truncated Newton method with a 245e76c431SBarry Smith line search. By default a cubic backtracking line search 255e76c431SBarry Smith is employed, as described in the text "Numerical Methods for 265e76c431SBarry Smith Unconstrained Optimization and Nonlinear Equations" by Dennis 275e42470aSBarry Smith and Schnabel. See the examples in src/snes/examples. 285e76c431SBarry Smith */ 295e76c431SBarry Smith 305e42470aSBarry Smith int SNESSolve_LS(SNES snes,int *outits) 315e76c431SBarry Smith { 325e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 33761aaf1bSLois Curfman McInnes int maxits, i, history_len, ierr, lits, lsfail; 34df60cc22SBarry Smith MatStructure flg = ALLMAT_DIFFERENT_NONZERO_PATTERN; 356b5873e3SBarry Smith double fnorm, gnorm, xnorm, ynorm, *history; 365e42470aSBarry Smith Vec Y, X, F, G, W, TMP; 375e76c431SBarry Smith 385e42470aSBarry Smith history = snes->conv_hist; /* convergence history */ 395e42470aSBarry Smith history_len = snes->conv_hist_len; /* convergence history length */ 405e42470aSBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 415e42470aSBarry Smith X = snes->vec_sol; /* solution vector */ 4239e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 435e42470aSBarry Smith Y = snes->work[0]; /* work vectors */ 445e42470aSBarry Smith G = snes->work[1]; 455e42470aSBarry Smith W = snes->work[2]; 465e76c431SBarry Smith 4778b31e54SBarry Smith ierr = SNESComputeInitialGuess(snes,X); CHKERRQ(ierr); /* X <- X_0 */ 485e42470aSBarry Smith VecNorm(X,&xnorm); /* xnorm = || X || */ 4978b31e54SBarry Smith ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr); /* (+/-) F(X) */ 505e42470aSBarry Smith VecNorm(F,&fnorm); /* fnorm <- ||F|| */ 515e42470aSBarry Smith snes->norm = fnorm; 525e76c431SBarry Smith if (history && history_len > 0) history[0] = fnorm; 53bbb6d6a8SBarry Smith if (snes->monitor) 54bbb6d6a8SBarry Smith {ierr = (*snes->monitor)(snes,0,fnorm,snes->monP); CHKERRQ(ierr);} 555e76c431SBarry Smith 565e76c431SBarry Smith for ( i=0; i<maxits; i++ ) { 575e42470aSBarry Smith snes->iter = i+1; 585e76c431SBarry Smith 595e76c431SBarry Smith /* Solve J Y = -F, where J is Jacobian matrix */ 60df60cc22SBarry Smith ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre, 61bbb6d6a8SBarry Smith &flg); CHKERRQ(ierr); 62a935fc98SLois Curfman McInnes ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre, 63a935fc98SLois Curfman McInnes flg); CHKERRQ(ierr); 6478b31e54SBarry Smith ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr); 6581b6cf68SLois Curfman McInnes ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr); 66761aaf1bSLois Curfman McInnes ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); 67761aaf1bSLois Curfman McInnes if (lsfail) snes->nfailures++; 6878b31e54SBarry Smith CHKERRQ(ierr); 695e76c431SBarry Smith 7039e2f89bSBarry Smith TMP = F; F = G; snes->vec_func_always = F; G = TMP; 7139e2f89bSBarry Smith TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 725e76c431SBarry Smith fnorm = gnorm; 735e76c431SBarry Smith 745e42470aSBarry Smith snes->norm = fnorm; 755e76c431SBarry Smith if (history && history_len > i+1) history[i+1] = fnorm; 765e42470aSBarry Smith VecNorm(X,&xnorm); /* xnorm = || X || */ 77bbb6d6a8SBarry Smith if (snes->monitor) 78bbb6d6a8SBarry Smith {(*snes->monitor)(snes,i+1,fnorm,snes->monP); CHKERRQ(ierr);} 795e76c431SBarry Smith 805e76c431SBarry Smith /* Test for convergence */ 81bbb6d6a8SBarry Smith if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) { 8239e2f89bSBarry Smith if (X != snes->vec_sol) { 8339e2f89bSBarry Smith VecCopy(X,snes->vec_sol); 8439e2f89bSBarry Smith snes->vec_sol_always = snes->vec_sol; 8539e2f89bSBarry Smith snes->vec_func_always = snes->vec_func; 8639e2f89bSBarry Smith } 875e76c431SBarry Smith break; 885e76c431SBarry Smith } 895e76c431SBarry Smith } 90*52392280SLois Curfman McInnes if (i == maxits) { 91*52392280SLois Curfman McInnes PLogInfo((PetscObject)snes, 92*52392280SLois Curfman McInnes "Maximum number of iterations has been reached: %d\n",maxits); 93*52392280SLois Curfman McInnes i--; 94*52392280SLois Curfman McInnes } 955e42470aSBarry Smith *outits = i+1; 965e42470aSBarry Smith return 0; 975e76c431SBarry Smith } 985e76c431SBarry Smith /* ------------------------------------------------------------ */ 995e42470aSBarry Smith int SNESSetUp_LS(SNES snes ) 1005e76c431SBarry Smith { 1015e42470aSBarry Smith int ierr; 10281b6cf68SLois Curfman McInnes snes->nwork = 4; 10378b31e54SBarry Smith ierr = VecGetVecs(snes->vec_sol,snes->nwork,&snes->work); CHKERRQ(ierr); 1045e42470aSBarry Smith PLogObjectParents(snes,snes->nwork,snes->work ); 10581b6cf68SLois Curfman McInnes snes->vec_sol_update_always = snes->work[3]; 1065e42470aSBarry Smith return 0; 1075e76c431SBarry Smith } 1085e76c431SBarry Smith /* ------------------------------------------------------------ */ 1095e42470aSBarry Smith int SNESDestroy_LS(PetscObject obj) 1105e76c431SBarry Smith { 1115e42470aSBarry Smith SNES snes = (SNES) obj; 1125e42470aSBarry Smith VecFreeVecs(snes->work, snes->nwork ); 11378b31e54SBarry Smith PETSCFREE(snes->data); 1145e42470aSBarry Smith return 0; 1155e76c431SBarry Smith } 1165e76c431SBarry Smith /* ------------------------------------------------------------ */ 1175e76c431SBarry Smith /*ARGSUSED*/ 1185e76c431SBarry Smith /*@ 1195e42470aSBarry Smith SNESNoLineSearch - This routine is not a line search at all; 1205e76c431SBarry Smith it simply uses the full Newton step. Thus, this routine is intended 1215e76c431SBarry Smith to serve as a template and is not recommended for general use. 1225e76c431SBarry Smith 1235e76c431SBarry Smith Input Parameters: 1245e42470aSBarry Smith . snes - nonlinear context 1255e76c431SBarry Smith . x - current iterate 1265e76c431SBarry Smith . f - residual evaluated at x 1275e76c431SBarry Smith . y - search direction (contains new iterate on output) 1285e76c431SBarry Smith . w - work vector 1295e76c431SBarry Smith . fnorm - 2-norm of f 1305e76c431SBarry Smith 131c4a48953SLois Curfman McInnes Output Parameters: 1325e76c431SBarry Smith . g - residual evaluated at new iterate y 1335e76c431SBarry Smith . y - new iterate (contains search direction on input) 1345e76c431SBarry Smith . gnorm - 2-norm of g 1355e76c431SBarry Smith . ynorm - 2-norm of search length 136761aaf1bSLois Curfman McInnes . flag - set to 0, indicating a successful line search 1375e76c431SBarry Smith 138c4a48953SLois Curfman McInnes Options Database Key: 139c4a48953SLois Curfman McInnes $ -snes_line_search basic 140c4a48953SLois Curfman McInnes 14128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 14228ae5a21SLois Curfman McInnes 143f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 144a935fc98SLois Curfman McInnes SNESSetLineSearchRoutine() 1455e76c431SBarry Smith @*/ 1465e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 147761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag ) 1485e76c431SBarry Smith { 1495e42470aSBarry Smith int ierr; 1505e42470aSBarry Smith Scalar one = 1.0; 151761aaf1bSLois Curfman McInnes *flag = 0; 1527857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 1535e42470aSBarry Smith VecNorm(y,ynorm); /* ynorm = || y || */ 1545e42470aSBarry Smith VecAXPY(&one,x,y); /* y <- x + y */ 15578b31e54SBarry Smith ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); 1565e42470aSBarry Smith VecNorm(g,gnorm); /* gnorm = || g || */ 1577857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 158761aaf1bSLois Curfman McInnes return 0; 1595e76c431SBarry Smith } 1605e76c431SBarry Smith /* ------------------------------------------------------------------ */ 1615e76c431SBarry Smith /*@ 162f59ffdedSLois Curfman McInnes SNESCubicLineSearch - This routine performs a cubic line search and 163f59ffdedSLois Curfman McInnes is the default line search method. 1645e76c431SBarry Smith 1655e76c431SBarry Smith Input Parameters: 1665e42470aSBarry Smith . snes - nonlinear context 1675e76c431SBarry Smith . x - current iterate 1685e76c431SBarry Smith . f - residual evaluated at x 1695e76c431SBarry Smith . y - search direction (contains new iterate on output) 1705e76c431SBarry Smith . w - work vector 1715e76c431SBarry Smith . fnorm - 2-norm of f 1725e76c431SBarry Smith 1735e76c431SBarry Smith Output parameters: 1745e76c431SBarry Smith . g - residual evaluated at new iterate y 1755e76c431SBarry Smith . y - new iterate (contains search direction on input) 1765e76c431SBarry Smith . gnorm - 2-norm of g 1775e76c431SBarry Smith . ynorm - 2-norm of search length 178761aaf1bSLois Curfman McInnes . flag - 0 if line search succeeds; -1 on failure. 1795e76c431SBarry Smith 180c4a48953SLois Curfman McInnes Options Database Key: 181c4a48953SLois Curfman McInnes $ -snes_line_search cubic 182c4a48953SLois Curfman McInnes 1835e76c431SBarry Smith Notes: 1845e76c431SBarry Smith This line search is taken from "Numerical Methods for Unconstrained 1855e76c431SBarry Smith Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 1865e76c431SBarry Smith 18728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic 18828ae5a21SLois Curfman McInnes 189f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine() 1905e76c431SBarry Smith @*/ 1915e42470aSBarry Smith int SNESCubicLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 192761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 1935e76c431SBarry Smith { 1945e42470aSBarry Smith double steptol, initslope; 1955e42470aSBarry Smith double lambdaprev, gnormprev; 1965e76c431SBarry Smith double a, b, d, t1, t2; 1976b5873e3SBarry Smith #if defined(PETSC_COMPLEX) 1985e42470aSBarry Smith Scalar cinitslope,clambda; 1996b5873e3SBarry Smith #endif 2005e42470aSBarry Smith int ierr,count; 2015e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 2025e42470aSBarry Smith Scalar one = 1.0,scale; 2035e42470aSBarry Smith double maxstep,minlambda,alpha,lambda,lambdatemp; 2045e76c431SBarry Smith 2057857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 206761aaf1bSLois Curfman McInnes *flag = 0; 2075e76c431SBarry Smith alpha = neP->alpha; 2085e76c431SBarry Smith maxstep = neP->maxstep; 2095e76c431SBarry Smith steptol = neP->steptol; 2105e76c431SBarry Smith 2115e42470aSBarry Smith VecNorm(y, ynorm ); 2125e76c431SBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 2135e42470aSBarry Smith scale = maxstep/(*ynorm); 2146b5873e3SBarry Smith #if defined(PETSC_COMPLEX) 2156b5873e3SBarry Smith PLogInfo((PetscObject)snes,"Scaling step by %g\n",real(scale)); 2166b5873e3SBarry Smith #else 2175e42470aSBarry Smith PLogInfo((PetscObject)snes,"Scaling step by %g\n",scale); 2186b5873e3SBarry Smith #endif 2195e42470aSBarry Smith VecScale(&scale, y ); 2205e76c431SBarry Smith *ynorm = maxstep; 2215e76c431SBarry Smith } 2225e76c431SBarry Smith minlambda = steptol/(*ynorm); 2235e42470aSBarry Smith #if defined(PETSC_COMPLEX) 2245e42470aSBarry Smith VecDot(f, y, &cinitslope ); 2255e42470aSBarry Smith initslope = real(cinitslope); 2265e42470aSBarry Smith #else 2275e42470aSBarry Smith VecDot(f, y, &initslope ); 2285e42470aSBarry Smith #endif 2295e76c431SBarry Smith if (initslope > 0.0) initslope = -initslope; 2305e76c431SBarry Smith if (initslope == 0.0) initslope = -1.0; 2315e76c431SBarry Smith 2325e42470aSBarry Smith VecCopy(y, w ); 2335e42470aSBarry Smith VecAXPY(&one, x, w ); 23478b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 2355e42470aSBarry Smith VecNorm(g, gnorm ); 2365e76c431SBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* Sufficient reduction */ 2375e42470aSBarry Smith VecCopy(w, y ); 2385e42470aSBarry Smith PLogInfo((PetscObject)snes,"Using full step\n"); 239d93a2b8dSBarry Smith goto theend; 2405e76c431SBarry Smith } 2415e76c431SBarry Smith 2425e76c431SBarry Smith /* Fit points with quadratic */ 2435e76c431SBarry Smith lambda = 1.0; count = 0; 2445e76c431SBarry Smith lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope)); 2455e76c431SBarry Smith lambdaprev = lambda; 2465e76c431SBarry Smith gnormprev = *gnorm; 2475e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 2485e76c431SBarry Smith lambda = .1*lambda; 2495e76c431SBarry Smith } else lambda = lambdatemp; 2505e42470aSBarry Smith VecCopy(x, w ); 2515e42470aSBarry Smith #if defined(PETSC_COMPLEX) 2525e42470aSBarry Smith clambda = lambda; VecAXPY(&clambda, y, w ); 2535e42470aSBarry Smith #else 2545e42470aSBarry Smith VecAXPY(&lambda, y, w ); 2555e42470aSBarry Smith #endif 25678b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 2575e42470aSBarry Smith VecNorm(g, gnorm ); 2585e76c431SBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* sufficient reduction */ 2595e42470aSBarry Smith VecCopy(w, y ); 2605e42470aSBarry Smith PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda); 261d93a2b8dSBarry Smith goto theend; 2625e76c431SBarry Smith } 2635e76c431SBarry Smith 2645e76c431SBarry Smith /* Fit points with cubic */ 2655e76c431SBarry Smith count = 1; 2665e76c431SBarry Smith while (1) { 2675e76c431SBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 2685e42470aSBarry Smith PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count); 2695e42470aSBarry Smith PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n", 2705e76c431SBarry Smith fnorm,*gnorm, *ynorm,lambda); 2715e42470aSBarry Smith VecCopy(w, y ); 272761aaf1bSLois Curfman McInnes *flag = -1; break; 2735e76c431SBarry Smith } 2745e76c431SBarry Smith t1 = *gnorm - fnorm - lambda*initslope; 2755e76c431SBarry Smith t2 = gnormprev - fnorm - lambdaprev*initslope; 2765e76c431SBarry Smith a = (t1/(lambda*lambda) - 2775e76c431SBarry Smith t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 2785e76c431SBarry Smith b = (-lambdaprev*t1/(lambda*lambda) + 2795e76c431SBarry Smith lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 2805e76c431SBarry Smith d = b*b - 3*a*initslope; 2815e76c431SBarry Smith if (d < 0.0) d = 0.0; 2825e76c431SBarry Smith if (a == 0.0) { 2835e76c431SBarry Smith lambdatemp = -initslope/(2.0*b); 2845e76c431SBarry Smith } else { 2855e76c431SBarry Smith lambdatemp = (-b + sqrt(d))/(3.0*a); 2865e76c431SBarry Smith } 2875e76c431SBarry Smith if (lambdatemp > .5*lambda) { 2885e76c431SBarry Smith lambdatemp = .5*lambda; 2895e76c431SBarry Smith } 2905e76c431SBarry Smith lambdaprev = lambda; 2915e76c431SBarry Smith gnormprev = *gnorm; 2925e76c431SBarry Smith if (lambdatemp <= .1*lambda) { 2935e76c431SBarry Smith lambda = .1*lambda; 2945e76c431SBarry Smith } 2955e76c431SBarry Smith else lambda = lambdatemp; 2965e42470aSBarry Smith VecCopy( x, w ); 2975e42470aSBarry Smith #if defined(PETSC_COMPLEX) 2985e42470aSBarry Smith VecAXPY(&clambda, y, w ); 2995e42470aSBarry Smith #else 3005e42470aSBarry Smith VecAXPY(&lambda, y, w ); 3015e42470aSBarry Smith #endif 30278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 3035e42470aSBarry Smith VecNorm(g, gnorm ); 3045e76c431SBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* is reduction enough */ 3055e42470aSBarry Smith VecCopy(w, y ); 3065e42470aSBarry Smith PLogInfo((PetscObject)snes,"Cubically determined step, lambda %g\n",lambda); 307761aaf1bSLois Curfman McInnes *flag = -1; break; 3085e76c431SBarry Smith } 3095e76c431SBarry Smith count++; 3105e76c431SBarry Smith } 311d93a2b8dSBarry Smith theend: 3127857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 3135e42470aSBarry Smith return 0; 3145e76c431SBarry Smith } 315*52392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 3165e42470aSBarry Smith /*@ 3175e42470aSBarry Smith SNESQuadraticLineSearch - This routine performs a cubic line search. 3185e76c431SBarry Smith 3195e42470aSBarry Smith Input Parameters: 320f59ffdedSLois Curfman McInnes . snes - the SNES context 3215e42470aSBarry Smith . x - current iterate 3225e42470aSBarry Smith . f - residual evaluated at x 3235e42470aSBarry Smith . y - search direction (contains new iterate on output) 3245e42470aSBarry Smith . w - work vector 3255e42470aSBarry Smith . fnorm - 2-norm of f 3265e42470aSBarry Smith 327c4a48953SLois Curfman McInnes Output Parameters: 3285e42470aSBarry Smith . g - residual evaluated at new iterate y 3295e42470aSBarry Smith . y - new iterate (contains search direction on input) 3305e42470aSBarry Smith . gnorm - 2-norm of g 3315e42470aSBarry Smith . ynorm - 2-norm of search length 332761aaf1bSLois Curfman McInnes . flag - 0 if line search succeeds; -1 on failure. 3335e42470aSBarry Smith 334c4a48953SLois Curfman McInnes Options Database Key: 335c4a48953SLois Curfman McInnes $ -snes_line_search quadratic 336c4a48953SLois Curfman McInnes 3375e42470aSBarry Smith Notes: 338f59ffdedSLois Curfman McInnes Use SNESSetLineSearchRoutine() 339f59ffdedSLois Curfman McInnes to set this routine within the SNES_NLS method. 3405e42470aSBarry Smith 341f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search 342f59ffdedSLois Curfman McInnes 343f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine() 3445e42470aSBarry Smith @*/ 3455e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w, 346761aaf1bSLois Curfman McInnes double fnorm, double *ynorm, double *gnorm,int *flag) 3475e76c431SBarry Smith { 3485e42470aSBarry Smith double steptol, initslope; 3495e42470aSBarry Smith double lambdaprev, gnormprev; 3506b5873e3SBarry Smith #if defined(PETSC_COMPLEX) 3515e42470aSBarry Smith Scalar cinitslope,clambda; 3526b5873e3SBarry Smith #endif 3535e42470aSBarry Smith int ierr,count; 3545e42470aSBarry Smith SNES_LS *neP = (SNES_LS *) snes->data; 3555e42470aSBarry Smith Scalar one = 1.0,scale; 3565e42470aSBarry Smith double maxstep,minlambda,alpha,lambda,lambdatemp; 3575e76c431SBarry Smith 3587857610eSBarry Smith PLogEventBegin(SNES_LineSearch,snes,x,f,g); 359761aaf1bSLois Curfman McInnes *flag = 0; 3605e42470aSBarry Smith alpha = neP->alpha; 3615e42470aSBarry Smith maxstep = neP->maxstep; 3625e42470aSBarry Smith steptol = neP->steptol; 3635e76c431SBarry Smith 3645e42470aSBarry Smith VecNorm(y, ynorm ); 3655e42470aSBarry Smith if (*ynorm > maxstep) { /* Step too big, so scale back */ 3665e42470aSBarry Smith scale = maxstep/(*ynorm); 3675e42470aSBarry Smith VecScale(&scale, y ); 3685e42470aSBarry Smith *ynorm = maxstep; 3695e76c431SBarry Smith } 3705e42470aSBarry Smith minlambda = steptol/(*ynorm); 3715e42470aSBarry Smith #if defined(PETSC_COMPLEX) 3725e42470aSBarry Smith VecDot(f, y, &cinitslope ); 3735e42470aSBarry Smith initslope = real(cinitslope); 3745e42470aSBarry Smith #else 3755e42470aSBarry Smith VecDot(f, y, &initslope ); 3765e42470aSBarry Smith #endif 3775e42470aSBarry Smith if (initslope > 0.0) initslope = -initslope; 3785e42470aSBarry Smith if (initslope == 0.0) initslope = -1.0; 3795e42470aSBarry Smith 3805e42470aSBarry Smith VecCopy(y, w ); 3815e42470aSBarry Smith VecAXPY(&one, x, w ); 38278b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 3835e42470aSBarry Smith VecNorm(g, gnorm ); 3845e42470aSBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* Sufficient reduction */ 3855e42470aSBarry Smith VecCopy(w, y ); 3865e42470aSBarry Smith PLogInfo((PetscObject)snes,"Using full step\n"); 387d93a2b8dSBarry Smith goto theend; 3885e42470aSBarry Smith } 3895e42470aSBarry Smith 3905e42470aSBarry Smith /* Fit points with quadratic */ 3915e42470aSBarry Smith lambda = 1.0; count = 0; 3925e42470aSBarry Smith count = 1; 3935e42470aSBarry Smith while (1) { 3945e42470aSBarry Smith if (lambda <= minlambda) { /* bad luck; use full step */ 3955e42470aSBarry Smith PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count); 3965e42470aSBarry Smith PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n", 3975e42470aSBarry Smith fnorm,*gnorm, *ynorm,lambda); 3985e42470aSBarry Smith VecCopy(w, y ); 399761aaf1bSLois Curfman McInnes *flag = -1; break; 4005e42470aSBarry Smith } 4015e42470aSBarry Smith lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope)); 4025e42470aSBarry Smith lambdaprev = lambda; 4035e42470aSBarry Smith gnormprev = *gnorm; 4045e42470aSBarry Smith if (lambdatemp <= .1*lambda) { 4055e42470aSBarry Smith lambda = .1*lambda; 4065e42470aSBarry Smith } else lambda = lambdatemp; 4075e42470aSBarry Smith VecCopy(x, w ); 4085e42470aSBarry Smith #if defined(PETSC_COMPLEX) 4095e42470aSBarry Smith clambda = lambda; VecAXPY(&clambda, y, w ); 4105e42470aSBarry Smith #else 4115e42470aSBarry Smith VecAXPY(&lambda, y, w ); 4125e42470aSBarry Smith #endif 41378b31e54SBarry Smith ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr); 4145e42470aSBarry Smith VecNorm(g, gnorm ); 4155e42470aSBarry Smith if (*gnorm <= fnorm + alpha*initslope) { /* sufficient reduction */ 4165e42470aSBarry Smith VecCopy(w, y ); 4175e42470aSBarry Smith PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda); 41806259719SBarry Smith break; 4195e42470aSBarry Smith } 4205e42470aSBarry Smith count++; 4215e42470aSBarry Smith } 422d93a2b8dSBarry Smith theend: 4237857610eSBarry Smith PLogEventEnd(SNES_LineSearch,snes,x,f,g); 4245e42470aSBarry Smith return 0; 4255e76c431SBarry Smith } 4265e76c431SBarry Smith /* ------------------------------------------------------------ */ 427b1f0a012SBarry Smith /*@ 4285e42470aSBarry Smith SNESSetLineSearchRoutine - Sets the line search routine to be used 429fbe28522SBarry Smith by the method SNES_LS. 4305e76c431SBarry Smith 4315e76c431SBarry Smith Input Parameters: 432eafb4bcbSBarry Smith . snes - nonlinear context obtained from SNESCreate() 4335e76c431SBarry Smith . func - pointer to int function 4345e76c431SBarry Smith 435c4a48953SLois Curfman McInnes Available Routines: 436f59ffdedSLois Curfman McInnes . SNESCubicLineSearch() - default line search 437f59ffdedSLois Curfman McInnes . SNESQuadraticLineSearch() - quadratic line search 438f59ffdedSLois Curfman McInnes . SNESNoLineSearch() - the full Newton step (actually not a line search) 4395e76c431SBarry Smith 440c4a48953SLois Curfman McInnes Options Database Keys: 441c4a48953SLois Curfman McInnes $ -snes_line_search [basic,quadratic,cubic] 442c4a48953SLois Curfman McInnes 4435e76c431SBarry Smith Calling sequence of func: 444f59ffdedSLois Curfman McInnes func (SNES snes, Vec x, Vec f, Vec g, Vec y, 445761aaf1bSLois Curfman McInnes Vec w, double fnorm, double *ynorm, 446761aaf1bSLois Curfman McInnes double *gnorm, *flag) 4475e76c431SBarry Smith 4485e76c431SBarry Smith Input parameters for func: 4495e42470aSBarry Smith . snes - nonlinear context 4505e76c431SBarry Smith . x - current iterate 4515e76c431SBarry Smith . f - residual evaluated at x 4525e76c431SBarry Smith . y - search direction (contains new iterate on output) 4535e76c431SBarry Smith . w - work vector 4545e76c431SBarry Smith . fnorm - 2-norm of f 4555e76c431SBarry Smith 4565e76c431SBarry Smith Output parameters for func: 4575e76c431SBarry Smith . g - residual evaluated at new iterate y 4585e76c431SBarry Smith . y - new iterate (contains search direction on input) 4595e76c431SBarry Smith . gnorm - 2-norm of g 4605e76c431SBarry Smith . ynorm - 2-norm of search length 461761aaf1bSLois Curfman McInnes . flag - set to 0 if the line search succeeds; a nonzero integer 462761aaf1bSLois Curfman McInnes on failure. 463f59ffdedSLois Curfman McInnes 464f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine 465f59ffdedSLois Curfman McInnes 466f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch() 4675e76c431SBarry Smith @*/ 4685e42470aSBarry Smith int SNESSetLineSearchRoutine(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec, 469761aaf1bSLois Curfman McInnes double,double *,double*,int*) ) 4705e76c431SBarry Smith { 471fbe28522SBarry Smith if ((snes)->type == SNES_NLS) 4725e42470aSBarry Smith ((SNES_LS *)(snes->data))->LineSearch = func; 4735e42470aSBarry Smith return 0; 4745e76c431SBarry Smith } 475*52392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 4765e42470aSBarry Smith static int SNESPrintHelp_LS(SNES snes) 4775e42470aSBarry Smith { 4785e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 479*52392280SLois Curfman McInnes char *p; 480*52392280SLois Curfman McInnes if (snes->prefix) p = snes->prefix; else p = "-"; 481*52392280SLois Curfman McInnes MPIU_printf(snes->comm," method ls (system of nonlinear equations):\n"); 482*52392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search [basic,quadratic,cubic]\n",p); 483*52392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search_alpha alpha (default %g)\n",p,ls->alpha); 484*52392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search_maxstep max (default %g)\n",p,ls->maxstep); 485*52392280SLois Curfman McInnes MPIU_printf(snes->comm," %ssnes_line_search_steptol tol (default %g)\n",p,ls->steptol); 4866b5873e3SBarry Smith return 0; 4875e42470aSBarry Smith } 488*52392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 489a935fc98SLois Curfman McInnes static int SNESView_LS(PetscObject obj,Viewer viewer) 490a935fc98SLois Curfman McInnes { 491a935fc98SLois Curfman McInnes SNES snes = (SNES)obj; 492a935fc98SLois Curfman McInnes SNES_LS *ls = (SNES_LS *)snes->data; 493a935fc98SLois Curfman McInnes FILE *fd = ViewerFileGetPointer_Private(viewer); 494a935fc98SLois Curfman McInnes char *cstring; 495a935fc98SLois Curfman McInnes 496a935fc98SLois Curfman McInnes if (ls->LineSearch == SNESNoLineSearch) 497a935fc98SLois Curfman McInnes cstring = "SNESNoLineSearch"; 498a935fc98SLois Curfman McInnes else if (ls->LineSearch == SNESQuadraticLineSearch) 499a935fc98SLois Curfman McInnes cstring = "SNESQuadraticLineSearch"; 500a935fc98SLois Curfman McInnes else if (ls->LineSearch == SNESCubicLineSearch) 501a935fc98SLois Curfman McInnes cstring = "SNESCubicLineSearch"; 502a935fc98SLois Curfman McInnes else 503a935fc98SLois Curfman McInnes cstring = "unknown"; 504a935fc98SLois Curfman McInnes MPIU_fprintf(snes->comm,fd," line search variant: %s\n",cstring); 505a935fc98SLois Curfman McInnes MPIU_fprintf(snes->comm,fd," alpha=%g, maxstep=%g, steptol=%g\n", 506a935fc98SLois Curfman McInnes ls->alpha,ls->maxstep,ls->steptol); 507a935fc98SLois Curfman McInnes return 0; 508a935fc98SLois Curfman McInnes } 509*52392280SLois Curfman McInnes /* ------------------------------------------------------------------ */ 5105e42470aSBarry Smith static int SNESSetFromOptions_LS(SNES snes) 5115e42470aSBarry Smith { 5125e42470aSBarry Smith SNES_LS *ls = (SNES_LS *)snes->data; 5135e42470aSBarry Smith char ver[16]; 5145e42470aSBarry Smith double tmp; 5155e42470aSBarry Smith 516df60cc22SBarry Smith if (OptionsGetDouble(snes->prefix,"-snes_line_search_alpa",&tmp)) { 5175e42470aSBarry Smith ls->alpha = tmp; 5185e42470aSBarry Smith } 519df60cc22SBarry Smith if (OptionsGetDouble(snes->prefix,"-snes_line_search_maxstep",&tmp)) { 5205e42470aSBarry Smith ls->maxstep = tmp; 5215e42470aSBarry Smith } 522df60cc22SBarry Smith if (OptionsGetDouble(snes->prefix,"-snes_line_search_steptol",&tmp)) { 5235e42470aSBarry Smith ls->steptol = tmp; 5245e42470aSBarry Smith } 525df60cc22SBarry Smith if (OptionsGetString(snes->prefix,"-snes_line_search",ver,16)) { 5265e42470aSBarry Smith if (!strcmp(ver,"basic")) { 5275e42470aSBarry Smith SNESSetLineSearchRoutine(snes,SNESNoLineSearch); 5285e42470aSBarry Smith } 5295e42470aSBarry Smith else if (!strcmp(ver,"quadratic")) { 5305e42470aSBarry Smith SNESSetLineSearchRoutine(snes,SNESQuadraticLineSearch); 5315e42470aSBarry Smith } 5325e42470aSBarry Smith else if (!strcmp(ver,"cubic")) { 5335e42470aSBarry Smith SNESSetLineSearchRoutine(snes,SNESCubicLineSearch); 5345e42470aSBarry Smith } 5358c48e012SBarry Smith else {SETERRQ(1,"SNESSetFromOptions_LS: Unknown line search");} 5365e42470aSBarry Smith } 5375e42470aSBarry Smith return 0; 5385e42470aSBarry Smith } 5395e42470aSBarry Smith /* ------------------------------------------------------------ */ 5405e42470aSBarry Smith int SNESCreate_LS(SNES snes ) 5415e42470aSBarry Smith { 5425e42470aSBarry Smith SNES_LS *neP; 5435e42470aSBarry Smith 5441a3481a4SLois Curfman McInnes if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 5451a3481a4SLois Curfman McInnes "SNESCreate_LS: Valid for SNES_NONLINEAR_EQUATIONS problems only"); 54625c7317fSLois Curfman McInnes snes->type = SNES_EQ_NLS; 54749e3953dSBarry Smith snes->setup = SNESSetUp_LS; 54849e3953dSBarry Smith snes->solve = SNESSolve_LS; 5495e42470aSBarry Smith snes->destroy = SNESDestroy_LS; 550bbb6d6a8SBarry Smith snes->converged = SNESDefaultConverged; 55149e3953dSBarry Smith snes->printhelp = SNESPrintHelp_LS; 55249e3953dSBarry Smith snes->setfromoptions = SNESSetFromOptions_LS; 553a935fc98SLois Curfman McInnes snes->view = SNESView_LS; 5545e42470aSBarry Smith 55578b31e54SBarry Smith neP = PETSCNEW(SNES_LS); CHKPTRQ(neP); 5565e42470aSBarry Smith snes->data = (void *) neP; 5575e42470aSBarry Smith neP->alpha = 1.e-4; 5585e42470aSBarry Smith neP->maxstep = 1.e8; 5595e42470aSBarry Smith neP->steptol = 1.e-12; 5605e42470aSBarry Smith neP->LineSearch = SNESCubicLineSearch; 5615e42470aSBarry Smith return 0; 5625e42470aSBarry Smith } 5635e42470aSBarry Smith 5645e42470aSBarry Smith 5655e42470aSBarry Smith 5665e42470aSBarry Smith 567