$extrastylesheet
nonlinear_implicit_system.h
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 
00020 #ifndef LIBMESH_NONLINEAR_IMPLICIT_SYSTEM_H
00021 #define LIBMESH_NONLINEAR_IMPLICIT_SYSTEM_H
00022 
00023 // Local Includes
00024 #include "libmesh/implicit_system.h"
00025 
00026 // C++ includes
00027 
00028 namespace libMesh
00029 {
00030 
00031 
00032 // Forward declarations
00033 class DiffSolver;
00034 template<typename T> class NonlinearSolver;
00035 
00036 
00045 // ------------------------------------------------------------
00046 // NonlinearImplicitSystem class definition
00047 
00048 class NonlinearImplicitSystem : public ImplicitSystem
00049 {
00050 public:
00051 
00056   NonlinearImplicitSystem (EquationSystems& es,
00057                            const std::string& name,
00058                            const unsigned int number);
00059 
00063   virtual ~NonlinearImplicitSystem ();
00064 
00068   typedef NonlinearImplicitSystem sys_type;
00069 
00073   typedef ImplicitSystem Parent;
00074 
00079   class ComputeResidual
00080   {
00081   public:
00082     virtual ~ComputeResidual () {}
00087     virtual void residual (const NumericVector<Number>& X,
00088                            NumericVector<Number>& R,
00089                            sys_type& S) = 0;
00090   };
00091 
00092 
00097   class ComputeJacobian
00098   {
00099   public:
00100     virtual ~ComputeJacobian () {}
00101 
00106     virtual void jacobian (const NumericVector<Number>& X,
00107                            SparseMatrix<Number>& J,
00108                            sys_type& S) = 0;
00109   };
00110 
00111 
00116   class ComputeBounds
00117   {
00118   public:
00119     virtual ~ComputeBounds () {}
00120 
00125     virtual void bounds (NumericVector<Number>& XL,
00126                          NumericVector<Number>& XU,
00127                          sys_type& S) = 0;
00128   };
00129 
00140   class ComputeVectorSubspace
00141   {
00142   public:
00143     virtual ~ComputeVectorSubspace () {}
00144 
00150     virtual void operator()(std::vector<NumericVector<Number>*>&sp, sys_type& s);
00151   };
00152 
00157   class ComputeResidualandJacobian
00158   {
00159   public:
00160     virtual ~ComputeResidualandJacobian () {}
00161 
00167     virtual void residual_and_jacobian (const NumericVector<Number>& X,
00168                                         NumericVector<Number>* R,
00169                                         SparseMatrix<Number>*  J,
00170                                         sys_type& S) = 0;
00171   };
00172 
00176   sys_type & system () { return *this; }
00177 
00182   virtual void clear ();
00183 
00188   virtual void reinit ();
00189 
00193   virtual void solve ();
00194 
00200   virtual std::pair<unsigned int, Real>
00201   get_linear_solve_parameters() const;
00202 
00207   virtual void assembly(bool get_residual, bool get_jacobian,
00208                         bool apply_heterogeneous_constraints = false);
00209 
00214   virtual std::string system_type () const { return "NonlinearImplicit"; }
00215 
00222   UniquePtr<NonlinearSolver<Number> > nonlinear_solver;
00223 
00228   UniquePtr<DiffSolver> diff_solver;
00229 
00234   unsigned int n_nonlinear_iterations() const { return _n_nonlinear_iterations; }
00235 
00239   Real final_nonlinear_residual() const { return _final_nonlinear_residual; }
00240 
00246   unsigned get_current_nonlinear_iteration_number() const;
00247 
00248 
00249 protected:
00250 
00254   void set_solver_parameters();
00255 
00260   unsigned int _n_nonlinear_iterations;
00261 
00265   Real _final_nonlinear_residual;
00266 };
00267 
00268 
00269 
00270 } // namespace libMesh
00271 
00272 // ------------------------------------------------------------
00273 // NonlinearImplicitSystem inline methods
00274 
00275 
00276 #endif // LIBMESH_NONLINEAR_IMPLICIT_SYSTEM_H