$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 #ifndef LIBMESH_NONLINEAR_SOLVER_H 00021 #define LIBMESH_NONLINEAR_SOLVER_H 00022 00023 // Local includes 00024 #include "libmesh/libmesh_common.h" 00025 #include "libmesh/enum_solver_package.h" 00026 #include "libmesh/reference_counted_object.h" 00027 #include "libmesh/nonlinear_implicit_system.h" 00028 #include "libmesh/libmesh.h" 00029 #include "libmesh/parallel_object.h" 00030 #include "libmesh/auto_ptr.h" 00031 00032 // C++ includes 00033 #include <cstddef> 00034 00035 namespace libMesh 00036 { 00037 00038 // forward declarations 00039 template <typename T> class SparseMatrix; 00040 template <typename T> class NumericVector; 00041 template <typename T> class Preconditioner; 00042 00050 template <typename T> 00051 class NonlinearSolver : public ReferenceCountedObject<NonlinearSolver<T> >, 00052 public ParallelObject 00053 { 00054 public: 00058 typedef NonlinearImplicitSystem sys_type; 00059 00063 explicit 00064 NonlinearSolver (sys_type& s); 00065 00069 virtual ~NonlinearSolver (); 00070 00075 static UniquePtr<NonlinearSolver<T> > build(sys_type& s, 00076 const SolverPackage solver_package = libMesh::default_solver_package()); 00077 00082 bool initialized () const { return _is_initialized; } 00083 00087 virtual void clear () {} 00088 00093 virtual void init (const char* name = NULL) = 0; 00094 00098 virtual std::pair<unsigned int, Real> solve (SparseMatrix<T>&, // System Jacobian Matrix 00099 NumericVector<T>&, // Solution vector 00100 NumericVector<T>&, // Residual vector 00101 const double, // Stopping tolerance 00102 const unsigned int) = 0; // N. Iterations 00103 00108 virtual void print_converged_reason() { libmesh_not_implemented(); } 00109 00113 virtual int get_total_linear_iterations() = 0; 00114 00120 virtual unsigned get_current_nonlinear_iteration_number() const = 0; 00121 00126 void (* residual) (const NumericVector<Number>& X, 00127 NumericVector<Number>& R, 00128 sys_type& S); 00129 00134 NonlinearImplicitSystem::ComputeResidual *residual_object; 00135 00140 void (* jacobian) (const NumericVector<Number>& X, 00141 SparseMatrix<Number>& J, 00142 sys_type& S); 00143 00148 NonlinearImplicitSystem::ComputeJacobian *jacobian_object; 00149 00156 void (* matvec) (const NumericVector<Number>& X, 00157 NumericVector<Number>* R, 00158 SparseMatrix<Number>* J, 00159 sys_type& S); 00160 00167 NonlinearImplicitSystem::ComputeResidualandJacobian *residual_and_jacobian_object; 00168 00172 void (* bounds) (NumericVector<Number>& XL, 00173 NumericVector<Number>& XU, 00174 sys_type& S); 00178 NonlinearImplicitSystem::ComputeBounds *bounds_object; 00179 00186 void (* nullspace) (std::vector<NumericVector<Number>*>& sp, sys_type& S); 00187 00194 NonlinearImplicitSystem::ComputeVectorSubspace *nullspace_object; 00195 00201 void (* nearnullspace) (std::vector<NumericVector<Number>*>& sp, sys_type& S); 00202 00208 NonlinearImplicitSystem::ComputeVectorSubspace *nearnullspace_object; 00209 00210 00211 void (* user_presolve)(sys_type& S); 00212 00216 const sys_type & system () const { return _system; } 00217 00221 sys_type & system () { return _system; } 00222 00226 void attach_preconditioner(Preconditioner<T> * preconditioner); 00227 00231 unsigned int max_nonlinear_iterations; 00232 00236 unsigned int max_function_evaluations; 00237 00248 Real absolute_residual_tolerance; 00249 Real relative_residual_tolerance; 00250 00262 Real absolute_step_tolerance; 00263 Real relative_step_tolerance; 00264 00269 unsigned int max_linear_iterations; 00270 00275 Real initial_linear_tolerance; 00276 00280 Real minimum_linear_tolerance; 00281 00286 bool converged; 00287 00288 protected: 00292 sys_type& _system; 00293 00297 bool _is_initialized; 00298 00302 Preconditioner<T> * _preconditioner; 00303 }; 00304 00305 00306 00307 00308 /*----------------------- inline functions ----------------------------------*/ 00309 template <typename T> 00310 inline 00311 NonlinearSolver<T>::NonlinearSolver (sys_type& s) : 00312 ParallelObject (s), 00313 residual (NULL), 00314 residual_object (NULL), 00315 jacobian (NULL), 00316 jacobian_object (NULL), 00317 matvec (NULL), 00318 residual_and_jacobian_object (NULL), 00319 bounds (NULL), 00320 bounds_object (NULL), 00321 nullspace (NULL), 00322 nullspace_object (NULL), 00323 nearnullspace (NULL), 00324 nearnullspace_object (NULL), 00325 user_presolve (NULL), 00326 max_nonlinear_iterations(0), 00327 max_function_evaluations(0), 00328 absolute_residual_tolerance(0), 00329 relative_residual_tolerance(0), 00330 absolute_step_tolerance(0), 00331 relative_step_tolerance(0), 00332 max_linear_iterations(0), 00333 initial_linear_tolerance(0), 00334 minimum_linear_tolerance(0), 00335 converged(false), 00336 _system(s), 00337 _is_initialized (false), 00338 _preconditioner (NULL) 00339 { 00340 } 00341 00342 00343 00344 template <typename T> 00345 inline 00346 NonlinearSolver<T>::~NonlinearSolver () 00347 { 00348 this->clear (); 00349 } 00350 00351 00352 } // namespace libMesh 00353 00354 00355 #endif // LIBMESH_NONLINEAR_SOLVER_H