$extrastylesheet
00001 00002 // The libMesh Finite Element Library. 00003 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00004 00005 // This library is free software; you can redistribute it and/or 00006 // modify it under the terms of the GNU Lesser General Public 00007 // License as published by the Free Software Foundation; either 00008 // version 2.1 of the License, or (at your option) any later version. 00009 00010 // This library is distributed in the hope that it will be useful, 00011 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00012 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00013 // Lesser General Public License for more details. 00014 00015 // You should have received a copy of the GNU Lesser General Public 00016 // License along with this library; if not, write to the Free Software 00017 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00018 00019 00020 00021 #ifndef LIBMESH_DIFF_PHYSICS_H 00022 #define LIBMESH_DIFF_PHYSICS_H 00023 00024 // Local Includes 00025 #include "libmesh/libmesh.h" 00026 #include "libmesh/auto_ptr.h" 00027 #include "libmesh/diff_context.h" 00028 00029 // C++ includes 00030 #include <vector> 00031 00032 namespace libMesh 00033 { 00034 00035 // Forward Declarations 00036 class System; 00037 class TimeSolver; 00038 00039 template <typename T> class NumericVector; 00040 00049 // ------------------------------------------------------------ 00050 // DifferentiablePhysics class definition 00051 00052 class DifferentiablePhysics 00053 { 00054 public: 00055 00060 DifferentiablePhysics () : 00061 compute_internal_sides (false), 00062 _mesh_sys (NULL), 00063 _mesh_x_var (libMesh::invalid_uint), 00064 _mesh_y_var (libMesh::invalid_uint), 00065 _mesh_z_var (libMesh::invalid_uint) 00066 {} 00067 00071 virtual ~DifferentiablePhysics (); 00072 00076 virtual UniquePtr<DifferentiablePhysics> clone_physics() = 0; 00077 00081 virtual void clear_physics (); 00082 00086 virtual void init_physics (const System& sys); 00087 00101 virtual bool element_time_derivative (bool request_jacobian, 00102 DiffContext &) { 00103 return request_jacobian; 00104 } 00105 00119 virtual bool element_constraint (bool request_jacobian, 00120 DiffContext &) { 00121 return request_jacobian; 00122 } 00123 00130 bool compute_internal_sides; 00131 00148 virtual bool side_time_derivative (bool request_jacobian, 00149 DiffContext &) { 00150 return request_jacobian; 00151 } 00152 00168 virtual bool side_constraint (bool request_jacobian, 00169 DiffContext &) { 00170 return request_jacobian; 00171 } 00172 00186 virtual bool nonlocal_time_derivative (bool request_jacobian, 00187 DiffContext &) { 00188 return request_jacobian; 00189 } 00190 00204 virtual bool nonlocal_constraint (bool request_jacobian, 00205 DiffContext &) { 00206 return request_jacobian; 00207 } 00208 00209 00221 virtual void time_evolving (unsigned int var) { 00222 if (_time_evolving.size() <= var) 00223 _time_evolving.resize(var+1, false); 00224 _time_evolving[var] = true; 00225 } 00226 00234 bool is_time_evolving (unsigned int var) const { 00235 return _time_evolving[var]; 00236 } 00237 00246 virtual bool eulerian_residual (bool request_jacobian, 00247 DiffContext &) { 00248 return request_jacobian; 00249 } 00250 00265 virtual bool mass_residual (bool request_jacobian, 00266 DiffContext &) { 00267 return request_jacobian; 00268 } 00269 00282 virtual bool side_mass_residual (bool request_jacobian, 00283 DiffContext &) { 00284 return request_jacobian; 00285 } 00286 00303 virtual bool nonlocal_mass_residual (bool request_jacobian, 00304 DiffContext &c); 00305 00306 /* 00307 * Prepares the result of a build_context() call for use. 00308 * 00309 * Most FEMSystem-based problems will need to reimplement this in order to 00310 * call FE::get_*() as their particular physics requires. 00311 */ 00312 virtual void init_context(DiffContext &) {} 00313 00337 virtual void set_mesh_system(System* sys); 00338 00344 const System* get_mesh_system() const; 00345 00350 System* get_mesh_system(); 00351 00366 virtual void set_mesh_x_var(unsigned int var); 00367 00372 unsigned int get_mesh_x_var() const; 00373 00378 virtual void set_mesh_y_var(unsigned int var); 00379 00384 unsigned int get_mesh_y_var() const; 00385 00390 virtual void set_mesh_z_var(unsigned int var); 00391 00396 unsigned int get_mesh_z_var() const; 00397 00403 bool _eulerian_time_deriv (bool request_jacobian, 00404 DiffContext&); 00405 00406 00407 protected: 00408 00412 System *_mesh_sys; 00413 00417 unsigned int _mesh_x_var, _mesh_y_var, _mesh_z_var; 00418 00423 std::vector<bool> _time_evolving; 00424 }; 00425 00426 00427 00428 // ------------------------------------------------------------ 00429 // DifferentiablePhysics inline methods 00430 00431 00432 inline 00433 void DifferentiablePhysics::set_mesh_system(System* sys) 00434 { 00435 // For now we assume that we're doing fully coupled mesh motion 00436 // if (sys && sys != this) 00437 // libmesh_not_implemented(); 00438 00439 // For the foreseeable future we'll assume that we keep these 00440 // Systems in the same EquationSystems 00441 // libmesh_assert_equal_to (&this->get_equation_systems(), 00442 // &sys->get_equation_systems()); 00443 00444 // And for the immediate future this code may not even work 00445 libmesh_experimental(); 00446 00447 _mesh_sys = sys; 00448 } 00449 00450 00451 00452 inline 00453 void DifferentiablePhysics::set_mesh_x_var (unsigned int var) 00454 { 00455 _mesh_x_var = var; 00456 } 00457 00458 00459 00460 inline 00461 void DifferentiablePhysics::set_mesh_y_var (unsigned int var) 00462 { 00463 _mesh_y_var = var; 00464 } 00465 00466 00467 00468 inline 00469 void DifferentiablePhysics::set_mesh_z_var (unsigned int var) 00470 { 00471 _mesh_z_var = var; 00472 } 00473 00474 00475 00476 inline 00477 const System* DifferentiablePhysics::get_mesh_system() const 00478 { 00479 return _mesh_sys; 00480 } 00481 00482 inline 00483 System* DifferentiablePhysics::get_mesh_system() 00484 { 00485 return _mesh_sys; 00486 } 00487 00488 inline 00489 unsigned int DifferentiablePhysics::get_mesh_x_var() const 00490 { 00491 return _mesh_x_var; 00492 } 00493 00494 inline 00495 unsigned int DifferentiablePhysics::get_mesh_y_var() const 00496 { 00497 return _mesh_y_var; 00498 } 00499 00500 inline 00501 unsigned int DifferentiablePhysics::get_mesh_z_var() const 00502 { 00503 return _mesh_z_var; 00504 } 00505 00506 00507 00508 } // namespace libMesh 00509 00510 00511 #endif // LIBMESH_DIFF_PHYSICS_H