nonlinear_implicit_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_NONLINEAR_IMPLICIT_SYSTEM_H
21 #define LIBMESH_NONLINEAR_IMPLICIT_SYSTEM_H
22 
23 // Local Includes
25 
26 // C++ includes
27 
28 namespace libMesh
29 {
30 
31 
32 // Forward declarations
33 class DiffSolver;
34 template<typename T> class NonlinearSolver;
35 
36 
45 // ------------------------------------------------------------
46 // NonlinearImplicitSystem class definition
47 
49 {
50 public:
51 
57  const std::string& name,
58  const unsigned int number);
59 
63  virtual ~NonlinearImplicitSystem ();
64 
69 
74 
80  {
81  public:
82  virtual ~ComputeResidual () {}
87  virtual void residual (const NumericVector<Number>& X,
89  sys_type& S) = 0;
90  };
91 
92 
98  {
99  public:
100  virtual ~ComputeJacobian () {}
101 
106  virtual void jacobian (const NumericVector<Number>& X,
108  sys_type& S) = 0;
109  };
110 
111 
117  {
118  public:
119  virtual ~ComputeBounds () {}
120 
125  virtual void bounds (NumericVector<Number>& XL,
127  sys_type& S) = 0;
128  };
129 
141  {
142  public:
144 
150  virtual void operator()(std::vector<NumericVector<Number>*>&sp, sys_type& s);
151  };
152 
158  {
159  public:
161 
167  virtual void residual_and_jacobian (const NumericVector<Number>& X,
170  sys_type& S) = 0;
171  };
172 
176  sys_type & system () { return *this; }
177 
182  virtual void clear ();
183 
188  virtual void reinit ();
189 
193  virtual void solve ();
194 
200  virtual std::pair<unsigned int, Real>
202 
207  virtual void assembly(bool get_residual, bool get_jacobian);
208 
213  virtual std::string system_type () const { return "NonlinearImplicit"; }
214 
222 
228 
233  unsigned int n_nonlinear_iterations() const { return _n_nonlinear_iterations; }
234 
239 
246 
247 
248 protected:
249 
253  void set_solver_parameters();
254 
260 
265 };
266 
267 
268 
269 } // namespace libMesh
270 
271 // ------------------------------------------------------------
272 // NonlinearImplicitSystem inline methods
273 
274 
275 #endif // LIBMESH_NONLINEAR_IMPLICIT_SYSTEM_H

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:06 UTC

Hosted By:
SourceForge.net Logo