$extrastylesheet
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 00020 // C++ includes 00021 00022 // Local includes 00023 #include "libmesh/nonlinear_implicit_system.h" 00024 #include "libmesh/diff_solver.h" 00025 #include "libmesh/equation_systems.h" 00026 #include "libmesh/libmesh_logging.h" 00027 #include "libmesh/nonlinear_solver.h" 00028 #include "libmesh/sparse_matrix.h" 00029 00030 namespace libMesh 00031 { 00032 00033 // ------------------------------------------------------------ 00034 // NonlinearImplicitSystem implementation 00035 NonlinearImplicitSystem::NonlinearImplicitSystem (EquationSystems& es, 00036 const std::string& name_in, 00037 const unsigned int number_in) : 00038 00039 Parent (es, name_in, number_in), 00040 nonlinear_solver (NonlinearSolver<Number>::build(*this)), 00041 diff_solver (), 00042 _n_nonlinear_iterations (0), 00043 _final_nonlinear_residual (1.e20) 00044 { 00045 // Set default parameters 00046 // These were chosen to match the Petsc defaults 00047 es.parameters.set<Real> ("linear solver tolerance") = 1e-5; 00048 es.parameters.set<Real> ("linear solver minimum tolerance") = 1e-5; 00049 es.parameters.set<unsigned int>("linear solver maximum iterations") = 10000; 00050 00051 es.parameters.set<unsigned int>("nonlinear solver maximum iterations") = 50; 00052 es.parameters.set<unsigned int>("nonlinear solver maximum function evaluations") = 10000; 00053 00054 es.parameters.set<Real>("nonlinear solver absolute residual tolerance") = 1e-35; 00055 es.parameters.set<Real>("nonlinear solver relative residual tolerance") = 1e-8; 00056 es.parameters.set<Real>("nonlinear solver absolute step tolerance") = 1e-8; 00057 es.parameters.set<Real>("nonlinear solver relative step tolerance") = 1e-8; 00058 } 00059 00060 00061 00062 NonlinearImplicitSystem::~NonlinearImplicitSystem () 00063 { 00064 // Clear data 00065 this->clear(); 00066 } 00067 00068 00069 00070 void NonlinearImplicitSystem::clear () 00071 { 00072 // clear the nonlinear solver 00073 nonlinear_solver->clear(); 00074 00075 // clear the parent data 00076 Parent::clear(); 00077 } 00078 00079 00080 00081 void NonlinearImplicitSystem::reinit () 00082 { 00083 // re-initialize the nonlinear solver interface 00084 nonlinear_solver->clear(); 00085 00086 if (diff_solver.get()) 00087 diff_solver->reinit(); 00088 00089 // initialize parent data 00090 Parent::reinit(); 00091 } 00092 00093 00094 00095 void NonlinearImplicitSystem::set_solver_parameters () 00096 { 00097 // Get a reference to the EquationSystems 00098 const EquationSystems& es = 00099 this->get_equation_systems(); 00100 00101 // Get the user-specifiied nonlinear solver tolerances 00102 const unsigned int maxits = 00103 es.parameters.get<unsigned int>("nonlinear solver maximum iterations"); 00104 00105 const unsigned int maxfuncs = 00106 es.parameters.get<unsigned int>("nonlinear solver maximum function evaluations"); 00107 00108 const Real abs_resid_tol = 00109 es.parameters.get<Real>("nonlinear solver absolute residual tolerance"); 00110 00111 const Real rel_resid_tol = 00112 es.parameters.get<Real>("nonlinear solver relative residual tolerance"); 00113 00114 const Real abs_step_tol = 00115 es.parameters.get<Real>("nonlinear solver absolute step tolerance"); 00116 00117 const Real rel_step_tol = 00118 es.parameters.get<Real>("nonlinear solver relative step tolerance"); 00119 00120 // Get the user-specified linear solver toleranaces 00121 const unsigned int maxlinearits = 00122 es.parameters.get<unsigned int>("linear solver maximum iterations"); 00123 00124 const Real linear_tol = 00125 es.parameters.get<Real>("linear solver tolerance"); 00126 00127 const Real linear_min_tol = 00128 es.parameters.get<Real>("linear solver minimum tolerance"); 00129 00130 // Set all the parameters on the NonlinearSolver 00131 nonlinear_solver->max_nonlinear_iterations = maxits; 00132 nonlinear_solver->max_function_evaluations = maxfuncs; 00133 nonlinear_solver->absolute_residual_tolerance = abs_resid_tol; 00134 nonlinear_solver->relative_residual_tolerance = rel_resid_tol; 00135 nonlinear_solver->absolute_step_tolerance = abs_step_tol; 00136 nonlinear_solver->relative_step_tolerance = rel_step_tol; 00137 nonlinear_solver->max_linear_iterations = maxlinearits; 00138 nonlinear_solver->initial_linear_tolerance = linear_tol; 00139 nonlinear_solver->minimum_linear_tolerance = linear_min_tol; 00140 00141 if (diff_solver.get()) 00142 { 00143 diff_solver->max_nonlinear_iterations = maxits; 00144 diff_solver->absolute_residual_tolerance = abs_resid_tol; 00145 diff_solver->relative_residual_tolerance = rel_resid_tol; 00146 diff_solver->absolute_step_tolerance = abs_step_tol; 00147 diff_solver->relative_step_tolerance = rel_step_tol; 00148 diff_solver->max_linear_iterations = maxlinearits; 00149 diff_solver->initial_linear_tolerance = linear_tol; 00150 diff_solver->minimum_linear_tolerance = linear_min_tol; 00151 } 00152 } 00153 00154 00155 00156 void NonlinearImplicitSystem::solve () 00157 { 00158 // Log how long the nonlinear solve takes. 00159 START_LOG("solve()", "System"); 00160 00161 this->set_solver_parameters(); 00162 00163 if (diff_solver.get()) 00164 { 00165 diff_solver->solve(); 00166 00167 // Store the number of nonlinear iterations required to 00168 // solve and the final residual. 00169 _n_nonlinear_iterations = diff_solver->total_outer_iterations(); 00170 _final_nonlinear_residual = 0.; // FIXME - support this! 00171 } 00172 else 00173 { 00174 if (libMesh::on_command_line("--solver_system_names")) 00175 nonlinear_solver->init((this->name()+"_").c_str()); 00176 else 00177 nonlinear_solver->init(); 00178 00179 // Solve the nonlinear system. 00180 const std::pair<unsigned int, Real> rval = 00181 nonlinear_solver->solve (*matrix, *solution, *rhs, 00182 nonlinear_solver->relative_residual_tolerance, 00183 nonlinear_solver->max_linear_iterations); 00184 00185 // Store the number of nonlinear iterations required to 00186 // solve and the final residual. 00187 _n_nonlinear_iterations = rval.first; 00188 _final_nonlinear_residual = rval.second; 00189 } 00190 00191 // Stop logging the nonlinear solve 00192 STOP_LOG("solve()", "System"); 00193 00194 // Update the system after the solve 00195 this->update(); 00196 } 00197 00198 00199 00200 std::pair<unsigned int, Real> NonlinearImplicitSystem::get_linear_solve_parameters() const 00201 { 00202 if (diff_solver.get()) 00203 return std::make_pair(this->diff_solver->max_linear_iterations, 00204 this->diff_solver->relative_residual_tolerance); 00205 return std::make_pair(this->nonlinear_solver->max_linear_iterations, 00206 this->nonlinear_solver->relative_residual_tolerance); 00207 } 00208 00209 00210 00211 void NonlinearImplicitSystem::assembly(bool get_residual, 00212 bool get_jacobian, 00213 bool /*apply_heterogeneous_constraints*/) 00214 { 00215 // Get current_local_solution in sync 00216 this->update(); 00217 00218 //----------------------------------------------------------------------------- 00219 // if the user has provided both function pointers and objects only the pointer 00220 // will be used, so catch that as an error 00221 if (nonlinear_solver->jacobian && nonlinear_solver->jacobian_object) 00222 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Jacobian!"); 00223 00224 if (nonlinear_solver->residual && nonlinear_solver->residual_object) 00225 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Residual!"); 00226 00227 if (nonlinear_solver->matvec && nonlinear_solver->residual_and_jacobian_object) 00228 libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!"); 00229 00230 00231 if (get_jacobian) 00232 { 00233 if (nonlinear_solver->jacobian != NULL) 00234 nonlinear_solver->jacobian (*current_local_solution.get(), *matrix, *this); 00235 00236 else if (nonlinear_solver->jacobian_object != NULL) 00237 nonlinear_solver->jacobian_object->jacobian (*current_local_solution.get(), *matrix, *this); 00238 00239 else if (nonlinear_solver->matvec != NULL) 00240 nonlinear_solver->matvec (*current_local_solution.get(), get_residual?rhs:NULL, matrix, *this); 00241 00242 else if (nonlinear_solver->residual_and_jacobian_object != NULL) 00243 nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian (*current_local_solution.get(), get_residual?rhs:NULL, matrix, *this); 00244 00245 else 00246 libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!"); 00247 } 00248 00249 if (get_residual) 00250 { 00251 if (nonlinear_solver->residual != NULL) 00252 nonlinear_solver->residual (*current_local_solution.get(), *rhs, *this); 00253 00254 else if (nonlinear_solver->residual_object != NULL) 00255 nonlinear_solver->residual_object->residual (*current_local_solution.get(), *rhs, *this); 00256 00257 else if (nonlinear_solver->matvec != NULL) 00258 { 00259 // we might have already grabbed the residual and jacobian together 00260 if (!get_jacobian) 00261 nonlinear_solver->matvec (*current_local_solution.get(), rhs, NULL, *this); 00262 } 00263 00264 else if (nonlinear_solver->residual_and_jacobian_object != NULL) 00265 { 00266 // we might have already grabbed the residual and jacobian together 00267 if (!get_jacobian) 00268 nonlinear_solver->residual_and_jacobian_object->residual_and_jacobian (*current_local_solution.get(), rhs, NULL, *this); 00269 } 00270 00271 else 00272 libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!"); 00273 } 00274 else 00275 libmesh_assert(get_jacobian); // I can't believe you really wanted to assemble *nothing* 00276 } 00277 00278 00279 00280 00281 unsigned NonlinearImplicitSystem::get_current_nonlinear_iteration_number() const 00282 { 00283 return nonlinear_solver->get_current_nonlinear_iteration_number(); 00284 } 00285 00286 00287 00288 } // namespace libMesh