$extrastylesheet
newton_solver.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 #include "libmesh/diff_system.h"
00020 #include "libmesh/dof_map.h"
00021 #include "libmesh/libmesh_logging.h"
00022 #include "libmesh/linear_solver.h"
00023 #include "libmesh/newton_solver.h"
00024 #include "libmesh/numeric_vector.h"
00025 #include "libmesh/sparse_matrix.h"
00026 
00027 namespace libMesh
00028 {
00029 
00030 // SIGN from Numerical Recipes
00031 template <typename T>
00032 inline
00033 T SIGN(T a, T b)
00034 {
00035   return b >= 0 ? std::abs(a) : -std::abs(a);
00036 }
00037 
00038 Real NewtonSolver::line_search(Real tol,
00039                                Real last_residual,
00040                                Real &current_residual,
00041                                NumericVector<Number> &newton_iterate,
00042                                const NumericVector<Number> &linear_solution)
00043 {
00044   // Take a full step if we got a residual reduction or if we
00045   // aren't substepping
00046   if ((current_residual < last_residual) ||
00047       (!require_residual_reduction &&
00048        (!require_finite_residual || !libmesh_isnan(current_residual))))
00049     return 1.;
00050 
00051   // The residual vector
00052   NumericVector<Number> &rhs = *(_system.rhs);
00053 
00054   Real ax = 0.;  // First abscissa, don't take negative steps
00055   Real cx = 1.;  // Second abscissa, don't extrapolate steps
00056 
00057   // Find bx, a step length that gives lower residual than ax or cx
00058   Real bx = 1.;
00059 
00060   while (libmesh_isnan(current_residual) ||
00061          (current_residual > last_residual &&
00062           require_residual_reduction))
00063     {
00064       // Reduce step size to 1/2, 1/4, etc.
00065       Real substepdivision;
00066       if (brent_line_search && !libmesh_isnan(current_residual))
00067         {
00068           substepdivision = std::min(0.5, last_residual/current_residual);
00069           substepdivision = std::max(substepdivision, tol*2.);
00070         }
00071       else
00072         substepdivision = 0.5;
00073 
00074       newton_iterate.add (bx * (1.-substepdivision),
00075                           linear_solution);
00076       newton_iterate.close();
00077       bx *= substepdivision;
00078       if (verbose)
00079         libMesh::out << "  Shrinking Newton step to "
00080                      << bx << std::endl;
00081 
00082       // Check residual with fractional Newton step
00083       _system.assembly (true, false);
00084 
00085       rhs.close();
00086       current_residual = rhs.l2_norm();
00087       if (verbose)
00088         libMesh::out << "  Current Residual: "
00089                      << current_residual << std::endl;
00090 
00091       if (bx/2. < minsteplength &&
00092           (libmesh_isnan(current_residual) ||
00093            (current_residual > last_residual)))
00094         {
00095           libMesh::out << "Inexact Newton step FAILED at step "
00096                        << _outer_iterations << std::endl;
00097 
00098           if (!continue_after_backtrack_failure)
00099             {
00100               libmesh_convergence_failure();
00101             }
00102           else
00103             {
00104               libMesh::out << "Continuing anyway ..." << std::endl;
00105               _solve_result = DiffSolver::DIVERGED_BACKTRACKING_FAILURE;
00106               return bx;
00107             }
00108         }
00109     } // end while (current_residual > last_residual)
00110 
00111   // Now return that reduced-residual step, or  use Brent's method to
00112   // find a more optimal step.
00113 
00114   if (!brent_line_search)
00115     return bx;
00116 
00117   // Brent's method adapted from Numerical Recipes in C, ch. 10.2
00118   Real e = 0.;
00119 
00120   Real x = bx, w = bx, v = bx;
00121 
00122   // Residuals at bx
00123   Real fx = current_residual,
00124     fw = current_residual,
00125     fv = current_residual;
00126 
00127   // Max iterations for Brent's method loop
00128   const unsigned int max_i = 20;
00129 
00130   // for golden ratio steps
00131   const Real golden_ratio = 1.-(std::sqrt(5.)-1.)/2.;
00132 
00133   for (unsigned int i=1; i <= max_i; i++)
00134     {
00135       Real xm = (ax+cx)*0.5;
00136       Real tol1 = tol * std::abs(x) + tol*tol;
00137       Real tol2 = 2.0 * tol1;
00138 
00139       // Test if we're done
00140       if (std::abs(x-xm) <= (tol2 - 0.5 * (cx - ax)))
00141         return x;
00142 
00143       Real d;
00144 
00145       // Construct a parabolic fit
00146       if (std::abs(e) > tol1)
00147         {
00148           Real r = (x-w)*(fx-fv);
00149           Real q = (x-v)*(fx-fw);
00150           Real p = (x-v)*q-(x-w)*r;
00151           q = 2. * (q-r);
00152           if (q > 0.)
00153             p = -p;
00154           else
00155             q = std::abs(q);
00156           if (std::abs(p) >= std::abs(0.5*q*e) ||
00157               p <= q * (ax-x) ||
00158               p >= q * (cx-x))
00159             {
00160               // Take a golden section step
00161               e = x >= xm ? ax-x : cx-x;
00162               d = golden_ratio * e;
00163             }
00164           else
00165             {
00166               // Take a parabolic fit step
00167               d = p/q;
00168               if (x+d-ax < tol2 || cx-(x+d) < tol2)
00169                 d = SIGN(tol1, xm - x);
00170             }
00171         }
00172       else
00173         {
00174           // Take a golden section step
00175           e = x >= xm ? ax-x : cx-x;
00176           d = golden_ratio * e;
00177         }
00178 
00179       Real u = std::abs(d) >= tol1 ? x+d : x + SIGN(tol1,d);
00180 
00181       // Assemble the residual at the new steplength u
00182       newton_iterate.add (bx - u, linear_solution);
00183       newton_iterate.close();
00184       bx = u;
00185       if (verbose)
00186         libMesh::out << "  Shrinking Newton step to "
00187                      << bx << std::endl;
00188 
00189       _system.assembly (true, false);
00190 
00191       rhs.close();
00192       Real fu = current_residual = rhs.l2_norm();
00193       if (verbose)
00194         libMesh::out << "  Current Residual: "
00195                      << fu << std::endl;
00196 
00197       if (fu <= fx)
00198         {
00199           if (u >= x)
00200             ax = x;
00201           else
00202             cx = x;
00203           v = w;   w = x;   x = u;
00204           fv = fw; fw = fx; fx = fu;
00205         }
00206       else
00207         {
00208           if (u < x)
00209             ax = u;
00210           else
00211             cx = u;
00212           if (fu <= fw || w == x)
00213             {
00214               v = w;   w = u;
00215               fv = fw; fw = fu;
00216             }
00217           else if (fu <= fv || v == x || v == w)
00218             {
00219               v = u;
00220               fv = fu;
00221             }
00222         }
00223     }
00224 
00225   if (!quiet)
00226     libMesh::out << "Warning!  Too many iterations used in Brent line search!"
00227                  << std::endl;
00228   return bx;
00229 }
00230 
00231 
00232 NewtonSolver::NewtonSolver (sys_type& s)
00233   : Parent(s),
00234     require_residual_reduction(true),
00235     require_finite_residual(true),
00236     brent_line_search(true),
00237     track_linear_convergence(false),
00238     minsteplength(1e-5),
00239     linear_tolerance_multiplier(1e-3),
00240     linear_solver(LinearSolver<Number>::build(s.comm()))
00241 {
00242 }
00243 
00244 
00245 
00246 NewtonSolver::~NewtonSolver ()
00247 {
00248 }
00249 
00250 
00251 
00252 void NewtonSolver::init()
00253 {
00254   Parent::init();
00255 
00256   if (libMesh::on_command_line("--solver_system_names"))
00257     linear_solver->init((_system.name()+"_").c_str());
00258   else
00259     linear_solver->init();
00260 
00261   linear_solver->init_names(_system);
00262 }
00263 
00264 
00265 
00266 void NewtonSolver::reinit()
00267 {
00268   Parent::reinit();
00269 
00270   linear_solver->clear();
00271 
00272   linear_solver->init_names(_system);
00273 }
00274 
00275 
00276 
00277 unsigned int NewtonSolver::solve()
00278 {
00279   START_LOG("solve()", "NewtonSolver");
00280 
00281   // Reset any prior solve result
00282   _solve_result = INVALID_SOLVE_RESULT;
00283 
00284   NumericVector<Number> &newton_iterate = *(_system.solution);
00285 
00286   UniquePtr<NumericVector<Number> > linear_solution_ptr = newton_iterate.zero_clone();
00287   NumericVector<Number> &linear_solution = *linear_solution_ptr;
00288   NumericVector<Number> &rhs = *(_system.rhs);
00289 
00290   newton_iterate.close();
00291   linear_solution.close();
00292   rhs.close();
00293 
00294 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00295   _system.get_dof_map().enforce_constraints_exactly(_system);
00296 #endif
00297 
00298   SparseMatrix<Number> &matrix = *(_system.matrix);
00299 
00300   // Set starting linear tolerance
00301   Real current_linear_tolerance = initial_linear_tolerance;
00302 
00303   // Start counting our linear solver steps
00304   _inner_iterations = 0;
00305 
00306   // Now we begin the nonlinear loop
00307   for (_outer_iterations=0; _outer_iterations<max_nonlinear_iterations;
00308        ++_outer_iterations)
00309     {
00310       if (verbose)
00311         libMesh::out << "Assembling the System" << std::endl;
00312 
00313       _system.assembly(true, true);
00314       rhs.close();
00315       Real current_residual = rhs.l2_norm();
00316 
00317       if (libmesh_isnan(current_residual))
00318         {
00319           libMesh::out << "  Nonlinear solver DIVERGED at step "
00320                        << _outer_iterations
00321                        << " with residual Not-a-Number"
00322                        << std::endl;
00323           libmesh_convergence_failure();
00324           continue;
00325         }
00326 
00327       if (current_residual == 0)
00328         {
00329           if (verbose)
00330             libMesh::out << "Linear solve unnecessary; residual = 0"
00331                          << std::endl;
00332 
00333           // We're not doing a solve, but other code may reuse this
00334           // matrix.
00335           matrix.close();
00336 
00337           if (absolute_residual_tolerance > 0)
00338             _solve_result |= CONVERGED_ABSOLUTE_RESIDUAL;
00339           if (relative_residual_tolerance > 0)
00340             _solve_result |= CONVERGED_RELATIVE_RESIDUAL;
00341           if (absolute_step_tolerance > 0)
00342             _solve_result |= CONVERGED_ABSOLUTE_STEP;
00343           if (relative_step_tolerance > 0)
00344             _solve_result |= CONVERGED_RELATIVE_STEP;
00345 
00346           break;
00347         }
00348 
00349       // Prepare to take incomplete steps
00350       Real last_residual = current_residual;
00351 
00352       max_residual_norm = std::max (current_residual,
00353                                     max_residual_norm);
00354 
00355       // Compute the l2 norm of the whole solution
00356       Real norm_total = newton_iterate.l2_norm();
00357 
00358       max_solution_norm = std::max(max_solution_norm, norm_total);
00359 
00360       if (verbose)
00361         libMesh::out << "Nonlinear Residual: "
00362                      << current_residual << std::endl;
00363 
00364       // Make sure our linear tolerance is low enough
00365       current_linear_tolerance = std::min (current_linear_tolerance,
00366                                            current_residual * linear_tolerance_multiplier);
00367 
00368       // But don't let it be too small
00369       if (current_linear_tolerance < minimum_linear_tolerance)
00370         {
00371           current_linear_tolerance = minimum_linear_tolerance;
00372         }
00373 
00374       // At this point newton_iterate is the current guess, and
00375       // linear_solution is now about to become the NEGATIVE of the next
00376       // Newton step.
00377 
00378       // Our best initial guess for the linear_solution is zero!
00379       linear_solution.zero();
00380 
00381       if (verbose)
00382         libMesh::out << "Linear solve starting, tolerance "
00383                      << current_linear_tolerance << std::endl;
00384 
00385       // Solve the linear system.
00386       const std::pair<unsigned int, Real> rval =
00387         linear_solver->solve (matrix, _system.request_matrix("Preconditioner"),
00388                               linear_solution, rhs, current_linear_tolerance,
00389                               max_linear_iterations);
00390 
00391       if (track_linear_convergence)
00392         {
00393           LinearConvergenceReason linear_c_reason = linear_solver->get_converged_reason();
00394 
00395           // Check if something went wrong during the linear solve
00396           if (linear_c_reason < 0)
00397             {
00398               // The linear solver failed somehow
00399               _solve_result |= DiffSolver::DIVERGED_LINEAR_SOLVER_FAILURE;
00400               // Print a message
00401               libMesh::out << "Linear solver failed during Newton step, dropping out."
00402                            << std::endl;
00403               break;
00404             }
00405         }
00406 
00407       // We may need to localize a parallel solution
00408       _system.update ();
00409       // The linear solver may not have fit our constraints exactly
00410 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00411       _system.get_dof_map().enforce_constraints_exactly
00412         (_system, &linear_solution, /* homogeneous = */ true);
00413 #endif
00414 
00415       const unsigned int linear_steps = rval.first;
00416       libmesh_assert_less_equal (linear_steps, max_linear_iterations);
00417       _inner_iterations += linear_steps;
00418 
00419       const bool linear_solve_finished =
00420         !(linear_steps == max_linear_iterations);
00421 
00422       if (verbose)
00423         libMesh::out << "Linear solve finished, step " << linear_steps
00424                      << ", residual " << rval.second
00425                      << std::endl;
00426 
00427       // Compute the l2 norm of the nonlinear update
00428       Real norm_delta = linear_solution.l2_norm();
00429 
00430       if (verbose)
00431         libMesh::out << "Trying full Newton step" << std::endl;
00432       // Take a full Newton step
00433       newton_iterate.add (-1., linear_solution);
00434       newton_iterate.close();
00435 
00436       if (this->linear_solution_monitor.get())
00437         {
00438           // Compute the l2 norm of the whole solution
00439           norm_total = newton_iterate.l2_norm();
00440           rhs.close();
00441           (*this->linear_solution_monitor)(linear_solution, norm_delta,
00442                                            newton_iterate, norm_total,
00443                                            rhs, rhs.l2_norm(), _outer_iterations);
00444         }
00445 
00446       // Check residual with full Newton step, if that's useful for determining
00447       // whether to line search, whether to quit early, or whether to die after
00448       // hitting our max iteration count
00449       if (this->require_residual_reduction ||
00450           this->require_finite_residual ||
00451           _outer_iterations+1 < max_nonlinear_iterations ||
00452           !continue_after_max_iterations)
00453         {
00454           _system.assembly(true, false);
00455 
00456           rhs.close();
00457           current_residual = rhs.l2_norm();
00458           if (verbose)
00459             libMesh::out << "  Current Residual: "
00460                          << current_residual << std::endl;
00461 
00462           // don't fiddle around if we've already converged
00463           if (test_convergence(current_residual, norm_delta,
00464                                linear_solve_finished &&
00465                                current_residual <= last_residual))
00466             {
00467               if (!quiet)
00468                 print_convergence(_outer_iterations, current_residual,
00469                                   norm_delta, linear_solve_finished &&
00470                                   current_residual <= last_residual);
00471               _outer_iterations++;
00472               break; // out of _outer_iterations for loop
00473             }
00474         }
00475 
00476       // since we're not converged, backtrack if necessary
00477       Real steplength =
00478         this->line_search(std::sqrt(TOLERANCE),
00479                           last_residual, current_residual,
00480                           newton_iterate, linear_solution);
00481       norm_delta *= steplength;
00482 
00483       // Check to see if backtracking failed,
00484       // and break out of the nonlinear loop if so...
00485       if (_solve_result == DiffSolver::DIVERGED_BACKTRACKING_FAILURE)
00486         {
00487           _outer_iterations++;
00488           break; // out of _outer_iterations for loop
00489         }
00490 
00491       if (_outer_iterations + 1 >= max_nonlinear_iterations)
00492         {
00493           libMesh::out << "  Nonlinear solver reached maximum step "
00494                        << max_nonlinear_iterations << ", latest evaluated residual "
00495                        << current_residual << std::endl;
00496           if (continue_after_max_iterations)
00497             {
00498               _solve_result = DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS;
00499               libMesh::out << "  Continuing..." << std::endl;
00500             }
00501           else
00502             {
00503               libmesh_convergence_failure();
00504             }
00505           continue;
00506         }
00507 
00508       // Compute the l2 norm of the whole solution
00509       norm_total = newton_iterate.l2_norm();
00510 
00511       max_solution_norm = std::max(max_solution_norm, norm_total);
00512 
00513       // Print out information for the
00514       // nonlinear iterations.
00515       if (verbose)
00516         libMesh::out << "  Nonlinear step: |du|/|u| = "
00517                      << norm_delta / norm_total
00518                      << ", |du| = " << norm_delta
00519                      << std::endl;
00520 
00521       // Terminate the solution iteration if the difference between
00522       // this iteration and the last is sufficiently small.
00523       if (test_convergence(current_residual, norm_delta / steplength,
00524                            linear_solve_finished))
00525         {
00526           if (!quiet)
00527             print_convergence(_outer_iterations, current_residual,
00528                               norm_delta / steplength,
00529                               linear_solve_finished);
00530           _outer_iterations++;
00531           break; // out of _outer_iterations for loop
00532         }
00533     } // end nonlinear loop
00534 
00535   // The linear solver may not have fit our constraints exactly
00536 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00537   _system.get_dof_map().enforce_constraints_exactly(_system);
00538 #endif
00539 
00540   // We may need to localize a parallel solution
00541   _system.update ();
00542 
00543   STOP_LOG("solve()", "NewtonSolver");
00544 
00545   // Make sure we are returning something sensible as the
00546   // _solve_result, except in the edge case where we weren't really asked to
00547   // solve.
00548   libmesh_assert (_solve_result != DiffSolver::INVALID_SOLVE_RESULT ||
00549                   !max_nonlinear_iterations);
00550 
00551   return _solve_result;
00552 }
00553 
00554 
00555 
00556 bool NewtonSolver::test_convergence(Real current_residual,
00557                                     Real step_norm,
00558                                     bool linear_solve_finished)
00559 {
00560   // We haven't converged unless we pass a convergence test
00561   bool has_converged = false;
00562 
00563   // Is our absolute residual low enough?
00564   if (current_residual < absolute_residual_tolerance)
00565     {
00566       _solve_result |= CONVERGED_ABSOLUTE_RESIDUAL;
00567       has_converged = true;
00568     }
00569 
00570   // Is our relative residual low enough?
00571   if ((current_residual / max_residual_norm) <
00572       relative_residual_tolerance)
00573     {
00574       _solve_result |= CONVERGED_RELATIVE_RESIDUAL;
00575       has_converged = true;
00576     }
00577 
00578   // For incomplete linear solves, it's not safe to test step sizes
00579   if (!linear_solve_finished)
00580     {
00581       return has_converged;
00582     }
00583 
00584   // Is our absolute Newton step size small enough?
00585   if (step_norm < absolute_step_tolerance)
00586     {
00587       _solve_result |= CONVERGED_ABSOLUTE_STEP;
00588       has_converged = true;
00589     }
00590 
00591   // Is our relative Newton step size small enough?
00592   if (step_norm / max_solution_norm <
00593       relative_step_tolerance)
00594     {
00595       _solve_result |= CONVERGED_RELATIVE_STEP;
00596       has_converged = true;
00597     }
00598 
00599   return has_converged;
00600 }
00601 
00602 
00603 void NewtonSolver::print_convergence(unsigned int step_num,
00604                                      Real current_residual,
00605                                      Real step_norm,
00606                                      bool linear_solve_finished)
00607 {
00608   // Is our absolute residual low enough?
00609   if (current_residual < absolute_residual_tolerance)
00610     {
00611       libMesh::out << "  Nonlinear solver converged, step " << step_num
00612                    << ", residual " << current_residual
00613                    << std::endl;
00614     }
00615   else if (absolute_residual_tolerance)
00616     {
00617       if (verbose)
00618         libMesh::out << "  Nonlinear solver current_residual "
00619                      << current_residual << " > "
00620                      << (absolute_residual_tolerance) << std::endl;
00621     }
00622 
00623   // Is our relative residual low enough?
00624   if ((current_residual / max_residual_norm) <
00625       relative_residual_tolerance)
00626     {
00627       libMesh::out << "  Nonlinear solver converged, step " << step_num
00628                    << ", residual reduction "
00629                    << current_residual / max_residual_norm
00630                    << " < " << relative_residual_tolerance
00631                    << std::endl;
00632     }
00633   else if (relative_residual_tolerance)
00634     {
00635       if (verbose)
00636         libMesh::out << "  Nonlinear solver relative residual "
00637                      << (current_residual / max_residual_norm)
00638                      << " > " << relative_residual_tolerance
00639                      << std::endl;
00640     }
00641 
00642   // For incomplete linear solves, it's not safe to test step sizes
00643   if (!linear_solve_finished)
00644     return;
00645 
00646   // Is our absolute Newton step size small enough?
00647   if (step_norm < absolute_step_tolerance)
00648     {
00649       libMesh::out << "  Nonlinear solver converged, step " << step_num
00650                    << ", absolute step size "
00651                    << step_norm
00652                    << " < " << absolute_step_tolerance
00653                    << std::endl;
00654     }
00655   else if (absolute_step_tolerance)
00656     {
00657       if (verbose)
00658         libMesh::out << "  Nonlinear solver absolute step size "
00659                      << step_norm
00660                      << " > " << absolute_step_tolerance
00661                      << std::endl;
00662     }
00663 
00664   // Is our relative Newton step size small enough?
00665   if (step_norm / max_solution_norm <
00666       relative_step_tolerance)
00667     {
00668       libMesh::out << "  Nonlinear solver converged, step " << step_num
00669                    << ", relative step size "
00670                    << (step_norm / max_solution_norm)
00671                    << " < " << relative_step_tolerance
00672                    << std::endl;
00673     }
00674   else if (relative_step_tolerance)
00675     {
00676       if (verbose)
00677         libMesh::out << "  Nonlinear solver relative step size "
00678                      << (step_norm / max_solution_norm)
00679                      << " > " << relative_step_tolerance
00680                      << std::endl;
00681     }
00682 }
00683 
00684 } // namespace libMesh