newton_solver.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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 ¤t_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 minsteplength(1e-5), 00238 linear_tolerance_multiplier(1e-3), 00239 linear_solver(LinearSolver<Number>::build()) 00240 { 00241 } 00242 00243 00244 00245 NewtonSolver::~NewtonSolver () 00246 { 00247 } 00248 00249 00250 00251 void NewtonSolver::reinit() 00252 { 00253 Parent::reinit(); 00254 00255 linear_solver->clear(); 00256 } 00257 00258 00259 00260 unsigned int NewtonSolver::solve() 00261 { 00262 START_LOG("solve()", "NewtonSolver"); 00263 00264 // Reset any prior solve result 00265 _solve_result = INVALID_SOLVE_RESULT; 00266 00267 NumericVector<Number> &newton_iterate = *(_system.solution); 00268 00269 AutoPtr<NumericVector<Number> > linear_solution_ptr = newton_iterate.zero_clone(); 00270 NumericVector<Number> &linear_solution = *linear_solution_ptr; 00271 NumericVector<Number> &rhs = *(_system.rhs); 00272 00273 newton_iterate.close(); 00274 linear_solution.close(); 00275 rhs.close(); 00276 00277 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00278 _system.get_dof_map().enforce_constraints_exactly(_system); 00279 #endif 00280 00281 SparseMatrix<Number> &matrix = *(_system.matrix); 00282 00283 // Prepare to take incomplete steps 00284 Real last_residual=0.; 00285 00286 // Set starting linear tolerance 00287 Real current_linear_tolerance = initial_linear_tolerance; 00288 00289 // Start counting our linear solver steps 00290 _inner_iterations = 0; 00291 00292 // Now we begin the nonlinear loop 00293 for (_outer_iterations=0; _outer_iterations<max_nonlinear_iterations; 00294 ++_outer_iterations) 00295 { 00296 if (verbose) 00297 libMesh::out << "Assembling the System" << std::endl; 00298 00299 _system.assembly(true, true); 00300 rhs.close(); 00301 Real current_residual = rhs.l2_norm(); 00302 last_residual = current_residual; 00303 00304 if (libmesh_isnan(current_residual)) 00305 { 00306 libMesh::out << " Nonlinear solver DIVERGED at step " 00307 << _outer_iterations 00308 << " with residual Not-a-Number" 00309 << std::endl; 00310 libmesh_convergence_failure(); 00311 continue; 00312 } 00313 00314 max_residual_norm = std::max (current_residual, 00315 max_residual_norm); 00316 00317 // Compute the l2 norm of the whole solution 00318 Real norm_total = newton_iterate.l2_norm(); 00319 00320 max_solution_norm = std::max(max_solution_norm, norm_total); 00321 00322 if (verbose) 00323 libMesh::out << "Nonlinear Residual: " 00324 << current_residual << std::endl; 00325 00326 // Make sure our linear tolerance is low enough 00327 current_linear_tolerance = std::min (current_linear_tolerance, 00328 current_residual * linear_tolerance_multiplier); 00329 00330 // But don't let it be too small 00331 if (current_linear_tolerance < minimum_linear_tolerance) 00332 { 00333 current_linear_tolerance = minimum_linear_tolerance; 00334 } 00335 00336 // At this point newton_iterate is the current guess, and 00337 // linear_solution is now about to become the NEGATIVE of the next 00338 // Newton step. 00339 00340 // Our best initial guess for the linear_solution is zero! 00341 linear_solution.zero(); 00342 00343 if (verbose) 00344 libMesh::out << "Linear solve starting, tolerance " 00345 << current_linear_tolerance << std::endl; 00346 00347 // Solve the linear system. 00348 const std::pair<unsigned int, Real> rval = 00349 linear_solver->solve (matrix, _system.request_matrix("Preconditioner"), 00350 linear_solution, rhs, current_linear_tolerance, 00351 max_linear_iterations); 00352 00353 // We may need to localize a parallel solution 00354 _system.update (); 00355 // The linear solver may not have fit our constraints exactly 00356 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00357 _system.get_dof_map().enforce_constraints_exactly 00358 (_system, &linear_solution, /* homogeneous = */ true); 00359 #endif 00360 00361 const unsigned int linear_steps = rval.first; 00362 libmesh_assert_less_equal (linear_steps, max_linear_iterations); 00363 _inner_iterations += linear_steps; 00364 00365 const bool linear_solve_finished = 00366 !(linear_steps == max_linear_iterations); 00367 00368 if (verbose) 00369 libMesh::out << "Linear solve finished, step " << linear_steps 00370 << ", residual " << rval.second 00371 << std::endl; 00372 00373 // Compute the l2 norm of the nonlinear update 00374 Real norm_delta = linear_solution.l2_norm(); 00375 00376 if (verbose) 00377 libMesh::out << "Trying full Newton step" << std::endl; 00378 // Take a full Newton step 00379 newton_iterate.add (-1., linear_solution); 00380 newton_iterate.close(); 00381 00382 if (this->linear_solution_monitor.get()) { 00383 // Compute the l2 norm of the whole solution 00384 norm_total = newton_iterate.l2_norm(); 00385 (*this->linear_solution_monitor)(linear_solution, norm_delta, newton_iterate, norm_total); 00386 } 00387 00388 // Check residual with full Newton step, if that's useful for determining 00389 // whether to line search, whether to quit early, or whether to die after 00390 // hitting our max iteration count 00391 if (this->require_residual_reduction || 00392 this->require_finite_residual || 00393 _outer_iterations+1 < max_nonlinear_iterations || 00394 !continue_after_max_iterations) 00395 { 00396 _system.assembly(true, false); 00397 00398 rhs.close(); 00399 current_residual = rhs.l2_norm(); 00400 if (verbose) 00401 libMesh::out << " Current Residual: " 00402 << current_residual << std::endl; 00403 00404 // don't fiddle around if we've already converged 00405 if (test_convergence(current_residual, norm_delta, 00406 linear_solve_finished && 00407 current_residual <= last_residual)) 00408 { 00409 if (!quiet) 00410 print_convergence(_outer_iterations, current_residual, 00411 norm_delta, linear_solve_finished && 00412 current_residual <= last_residual); 00413 _outer_iterations++; 00414 break; // out of _outer_iterations for loop 00415 } 00416 } 00417 00418 // since we're not converged, backtrack if necessary 00419 Real steplength = 1; 00420 this->line_search(std::sqrt(TOLERANCE), 00421 last_residual, current_residual, 00422 newton_iterate, linear_solution); 00423 norm_delta *= steplength; 00424 00425 // Check to see if backtracking failed, 00426 // and break out of the nonlinear loop if so... 00427 if (_solve_result == DiffSolver::DIVERGED_BACKTRACKING_FAILURE) 00428 { 00429 _outer_iterations++; 00430 break; // out of _outer_iterations for loop 00431 } 00432 00433 if (_outer_iterations + 1 >= max_nonlinear_iterations) 00434 { 00435 libMesh::out << " Nonlinear solver reached maximum step " 00436 << max_nonlinear_iterations << ", latest evaluated residual " 00437 << current_residual << std::endl; 00438 if (continue_after_max_iterations) 00439 { 00440 _solve_result = DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS; 00441 libMesh::out << " Continuing..." << std::endl; 00442 } 00443 else 00444 { 00445 libmesh_convergence_failure(); 00446 } 00447 continue; 00448 } 00449 00450 // Compute the l2 norm of the whole solution 00451 norm_total = newton_iterate.l2_norm(); 00452 00453 max_solution_norm = std::max(max_solution_norm, norm_total); 00454 00455 // Print out information for the 00456 // nonlinear iterations. 00457 if (verbose) 00458 libMesh::out << " Nonlinear step: |du|/|u| = " 00459 << norm_delta / norm_total 00460 << ", |du| = " << norm_delta 00461 << std::endl; 00462 00463 // Terminate the solution iteration if the difference between 00464 // this iteration and the last is sufficiently small. 00465 if (test_convergence(current_residual, norm_delta / steplength, 00466 linear_solve_finished)) 00467 { 00468 if (!quiet) 00469 print_convergence(_outer_iterations, current_residual, 00470 norm_delta / steplength, 00471 linear_solve_finished); 00472 _outer_iterations++; 00473 break; // out of _outer_iterations for loop 00474 } 00475 } // end nonlinear loop 00476 00477 // The linear solver may not have fit our constraints exactly 00478 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00479 _system.get_dof_map().enforce_constraints_exactly(_system); 00480 #endif 00481 00482 // We may need to localize a parallel solution 00483 _system.update (); 00484 00485 STOP_LOG("solve()", "NewtonSolver"); 00486 00487 // Make sure we are returning something sensible as the 00488 // _solve_result, except in the edge case where we weren't really asked to 00489 // solve. 00490 libmesh_assert (_solve_result != DiffSolver::INVALID_SOLVE_RESULT || 00491 !max_nonlinear_iterations); 00492 00493 return _solve_result; 00494 } 00495 00496 00497 00498 bool NewtonSolver::test_convergence(Real current_residual, 00499 Real step_norm, 00500 bool linear_solve_finished) 00501 { 00502 // We haven't converged unless we pass a convergence test 00503 bool has_converged = false; 00504 00505 // Is our absolute residual low enough? 00506 if (current_residual < absolute_residual_tolerance) 00507 { 00508 _solve_result |= CONVERGED_ABSOLUTE_RESIDUAL; 00509 has_converged = true; 00510 } 00511 00512 // Is our relative residual low enough? 00513 if ((current_residual / max_residual_norm) < 00514 relative_residual_tolerance) 00515 { 00516 _solve_result |= CONVERGED_RELATIVE_RESIDUAL; 00517 has_converged = true; 00518 } 00519 00520 // For incomplete linear solves, it's not safe to test step sizes 00521 if (!linear_solve_finished) 00522 { 00523 return has_converged; 00524 } 00525 00526 // Is our absolute Newton step size small enough? 00527 if (step_norm < absolute_step_tolerance) 00528 { 00529 _solve_result |= CONVERGED_ABSOLUTE_STEP; 00530 has_converged = true; 00531 } 00532 00533 // Is our relative Newton step size small enough? 00534 if (step_norm / max_solution_norm < 00535 relative_step_tolerance) 00536 { 00537 _solve_result |= CONVERGED_RELATIVE_STEP; 00538 has_converged = true; 00539 } 00540 00541 return has_converged; 00542 } 00543 00544 00545 void NewtonSolver::print_convergence(unsigned int step_num, 00546 Real current_residual, 00547 Real step_norm, 00548 bool linear_solve_finished) 00549 { 00550 // Is our absolute residual low enough? 00551 if (current_residual < absolute_residual_tolerance) 00552 { 00553 libMesh::out << " Nonlinear solver converged, step " << step_num 00554 << ", residual " << current_residual 00555 << std::endl; 00556 } 00557 else if (absolute_residual_tolerance) 00558 { 00559 if (verbose) 00560 libMesh::out << " Nonlinear solver current_residual " 00561 << current_residual << " > " 00562 << (absolute_residual_tolerance) << std::endl; 00563 } 00564 00565 // Is our relative residual low enough? 00566 if ((current_residual / max_residual_norm) < 00567 relative_residual_tolerance) 00568 { 00569 libMesh::out << " Nonlinear solver converged, step " << step_num 00570 << ", residual reduction " 00571 << current_residual / max_residual_norm 00572 << " < " << relative_residual_tolerance 00573 << std::endl; 00574 } 00575 else if (relative_residual_tolerance) 00576 { 00577 if (verbose) 00578 libMesh::out << " Nonlinear solver relative residual " 00579 << (current_residual / max_residual_norm) 00580 << " > " << relative_residual_tolerance 00581 << std::endl; 00582 } 00583 00584 // For incomplete linear solves, it's not safe to test step sizes 00585 if (!linear_solve_finished) 00586 return; 00587 00588 // Is our absolute Newton step size small enough? 00589 if (step_norm < absolute_step_tolerance) 00590 { 00591 libMesh::out << " Nonlinear solver converged, step " << step_num 00592 << ", absolute step size " 00593 << step_norm 00594 << " < " << absolute_step_tolerance 00595 << std::endl; 00596 } 00597 else if (absolute_step_tolerance) 00598 { 00599 if (verbose) 00600 libMesh::out << " Nonlinear solver absolute step size " 00601 << step_norm 00602 << " > " << absolute_step_tolerance 00603 << std::endl; 00604 } 00605 00606 // Is our relative Newton step size small enough? 00607 if (step_norm / max_solution_norm < 00608 relative_step_tolerance) 00609 { 00610 libMesh::out << " Nonlinear solver converged, step " << step_num 00611 << ", relative step size " 00612 << (step_norm / max_solution_norm) 00613 << " < " << relative_step_tolerance 00614 << std::endl; 00615 } 00616 else if (relative_step_tolerance) 00617 { 00618 if (verbose) 00619 libMesh::out << " Nonlinear solver relative step size " 00620 << (step_norm / max_solution_norm) 00621 << " > " << relative_step_tolerance 00622 << std::endl; 00623 } 00624 } 00625 00626 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: