diff_system.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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) = 0; 00128 00134 virtual void solve (); 00135 00140 virtual std::pair<unsigned int, Real> adjoint_solve (const QoISet& qoi_indices = QoISet()); 00141 00145 virtual AutoPtr<DifferentiablePhysics> clone_physics() 00146 { libmesh_error(); 00147 // dummy 00148 return AutoPtr<DifferentiablePhysics>(this); } 00149 00153 virtual AutoPtr<DifferentiableQoI> clone() 00154 { libmesh_error(); 00155 // dummy 00156 return AutoPtr<DifferentiableQoI>(this); } 00157 00163 const DifferentiablePhysics* get_physics() const 00164 { return this->_diff_physics; } 00165 00170 DifferentiablePhysics* get_physics() 00171 { return this->_diff_physics; } 00172 00176 void attach_physics( DifferentiablePhysics* physics_in ) 00177 { this->_diff_physics = (physics_in->clone_physics()).release(); 00178 this->_diff_physics->init_physics(*this);} 00179 00184 const DifferentiableQoI* get_qoi() const 00185 { return this->diff_qoi; } 00186 00191 DifferentiableQoI* get_qoi() 00192 { return this->diff_qoi; } 00193 00197 void attach_qoi( DifferentiableQoI* qoi_in ) 00198 { this->diff_qoi = (qoi_in->clone()).release(); 00199 // User needs to resize qoi system qoi accordingly 00200 this->diff_qoi->init_qoi( this->qoi );} 00201 00206 AutoPtr<TimeSolver> time_solver; 00207 00211 void set_time_solver(AutoPtr<TimeSolver> _time_solver) 00212 { 00213 time_solver = _time_solver; 00214 } 00215 00219 TimeSolver& get_time_solver(); 00220 00224 const TimeSolver& get_time_solver() const; 00225 00230 Real deltat; 00231 00240 virtual AutoPtr<DiffContext> build_context(); 00241 00246 virtual void postprocess (){} 00247 00251 virtual void element_postprocess (DiffContext &) {} 00252 00257 virtual void side_postprocess (DiffContext &) {} 00258 00263 bool postprocess_sides; 00264 00269 bool print_solution_norms; 00270 00275 bool print_solutions; 00276 00280 bool print_residual_norms; 00281 00285 bool print_residuals; 00286 00290 bool print_jacobian_norms; 00291 00295 bool print_jacobians; 00296 00300 bool print_element_jacobians; 00301 00302 protected: 00303 00309 DifferentiablePhysics *_diff_physics; 00310 00316 DifferentiableQoI* diff_qoi; 00317 00322 virtual void init_data (); 00323 }; 00324 00325 // -------------------------------------------------------------- 00326 // DifferentiableSystem inline methods 00327 inline 00328 TimeSolver& DifferentiableSystem::get_time_solver() 00329 { 00330 libmesh_assert(time_solver.get()); 00331 libmesh_assert_equal_to (&(time_solver->system()), this); 00332 return *time_solver; 00333 } 00334 00335 inline 00336 const TimeSolver& DifferentiableSystem::get_time_solver() const 00337 { 00338 libmesh_assert(time_solver.get()); 00339 libmesh_assert_equal_to (&(time_solver->system()), this); 00340 return *time_solver; 00341 } 00342 00343 } // namespace libMesh 00344 00345 00346 #endif // LIBMESH_DIFF_SYSTEM_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: