$extrastylesheet
diff_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_DIFF_SYSTEM_H
00021 #define LIBMESH_DIFF_SYSTEM_H
00022 
00023 // Local Includes
00024 #include "libmesh/auto_ptr.h"
00025 #include "libmesh/diff_context.h"
00026 #include "libmesh/diff_physics.h"
00027 #include "libmesh/diff_qoi.h"
00028 #include "libmesh/implicit_system.h"
00029 #include "libmesh/time_solver.h"
00030 
00031 // C++ includes
00032 
00033 namespace libMesh
00034 {
00035 
00036 // Forward Declarations
00037 class TimeSolver;
00038 
00039 template <typename T> class NumericVector;
00040 
00053 // ------------------------------------------------------------
00054 // DifferentiableSystem class definition
00055 
00056 class DifferentiableSystem : public ImplicitSystem,
00057                              public virtual DifferentiablePhysics,
00058                              public virtual DifferentiableQoI
00059 {
00060 public:
00061 
00066   DifferentiableSystem (EquationSystems& es,
00067                         const std::string& name,
00068                         const unsigned int number);
00069 
00073   virtual ~DifferentiableSystem ();
00074 
00078   typedef DifferentiableSystem sys_type;
00079 
00083   typedef ImplicitSystem Parent;
00084 
00089   virtual void clear ();
00090 
00095   virtual void reinit ();
00096 
00101   virtual void assemble ();
00102 
00107   virtual LinearSolver<Number> *get_linear_solver() const;
00108 
00114   virtual std::pair<unsigned int, Real>
00115   get_linear_solve_parameters() const;
00116 
00121   virtual void release_linear_solver(LinearSolver<Number> *) const;
00122 
00127   virtual void assembly (bool get_residual, bool get_jacobian,
00128                          bool apply_heterogeneous_constraints = false) = 0;
00129 
00135   virtual void solve ();
00136 
00141   virtual std::pair<unsigned int, Real> adjoint_solve (const QoISet& qoi_indices = QoISet());
00142 
00146   virtual UniquePtr<DifferentiablePhysics> clone_physics()
00147   {
00148     libmesh_not_implemented();
00149     // dummy
00150     return UniquePtr<DifferentiablePhysics>(this);
00151   }
00152 
00156   virtual UniquePtr<DifferentiableQoI> clone()
00157   {
00158     libmesh_not_implemented();
00159     // dummy
00160     return UniquePtr<DifferentiableQoI>(this);
00161   }
00162 
00168   const DifferentiablePhysics* get_physics() const
00169   { return this->_diff_physics; }
00170 
00175   DifferentiablePhysics* get_physics()
00176   { return this->_diff_physics; }
00177 
00181   void attach_physics( DifferentiablePhysics* physics_in )
00182   { this->_diff_physics = (physics_in->clone_physics()).release();
00183     this->_diff_physics->init_physics(*this);}
00184 
00189   const DifferentiableQoI* get_qoi() const
00190   { return this->diff_qoi; }
00191 
00196   DifferentiableQoI* get_qoi()
00197   { return this->diff_qoi; }
00198 
00202   void attach_qoi( DifferentiableQoI* qoi_in )
00203   { this->diff_qoi = (qoi_in->clone()).release();
00204     // User needs to resize qoi system qoi accordingly
00205     this->diff_qoi->init_qoi( this->qoi );}
00206 
00211   UniquePtr<TimeSolver> time_solver;
00212 
00219   void set_time_solver(UniquePtr<TimeSolver> _time_solver)
00220   {
00221     time_solver.reset(_time_solver.release());
00222   }
00223 
00227   TimeSolver& get_time_solver();
00228 
00232   const TimeSolver& get_time_solver() const;
00233 
00238   Real deltat;
00239 
00248   virtual UniquePtr<DiffContext> build_context();
00249 
00254   virtual void postprocess (){}
00255 
00259   virtual void element_postprocess (DiffContext &) {}
00260 
00265   virtual void side_postprocess (DiffContext &) {}
00266 
00271   bool postprocess_sides;
00272 
00277   bool print_solution_norms;
00278 
00283   bool print_solutions;
00284 
00288   bool print_residual_norms;
00289 
00293   bool print_residuals;
00294 
00298   bool print_jacobian_norms;
00299 
00303   bool print_jacobians;
00304 
00308   bool print_element_solutions;
00309 
00313   bool print_element_residuals;
00314 
00318   bool print_element_jacobians;
00319 
00320 protected:
00321 
00327   DifferentiablePhysics *_diff_physics;
00328 
00334   DifferentiableQoI* diff_qoi;
00335 
00340   virtual void init_data ();
00341 };
00342 
00343 // --------------------------------------------------------------
00344 // DifferentiableSystem inline methods
00345 inline
00346 TimeSolver& DifferentiableSystem::get_time_solver()
00347 {
00348   libmesh_assert(time_solver.get());
00349   libmesh_assert_equal_to (&(time_solver->system()), this);
00350   return *time_solver;
00351 }
00352 
00353 inline
00354 const TimeSolver& DifferentiableSystem::get_time_solver() const
00355 {
00356   libmesh_assert(time_solver.get());
00357   libmesh_assert_equal_to (&(time_solver->system()), this);
00358   return *time_solver;
00359 }
00360 
00361 } // namespace libMesh
00362 
00363 
00364 #endif // LIBMESH_DIFF_SYSTEM_H