diff_physics.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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_DIFF_PHYSICS_H
21 #define LIBMESH_DIFF_PHYSICS_H
22 
23 // Local Includes
24 #include "libmesh/libmesh.h"
25 #include "libmesh/auto_ptr.h" // deprecated
26 
27 // C++ includes
28 #include <vector>
29 #include <memory>
30 
31 namespace libMesh
32 {
33 
34 // Forward Declarations
35 class System;
36 class DiffContext;
37 
76 {
77 public:
78 
84  compute_internal_sides (false),
85  _mesh_sys (nullptr),
89  {}
90 
94  virtual ~DifferentiablePhysics ();
95 
99  virtual std::unique_ptr<DifferentiablePhysics> clone_physics() = 0;
100 
104  virtual void clear_physics ();
105 
109  virtual void init_physics (const System & sys);
110 
124  virtual bool element_time_derivative (bool request_jacobian,
125  DiffContext &) {
126  return request_jacobian;
127  }
128 
142  virtual bool element_constraint (bool request_jacobian,
143  DiffContext &) {
144  return request_jacobian;
145  }
146 
154 
171  virtual bool side_time_derivative (bool request_jacobian,
172  DiffContext &) {
173  return request_jacobian;
174  }
175 
191  virtual bool side_constraint (bool request_jacobian,
192  DiffContext &) {
193  return request_jacobian;
194  }
195 
209  virtual bool nonlocal_time_derivative (bool request_jacobian,
210  DiffContext &) {
211  return request_jacobian;
212  }
213 
227  virtual bool nonlocal_constraint (bool request_jacobian,
228  DiffContext &) {
229  return request_jacobian;
230  }
231 
232 
248 #ifdef LIBMESH_ENABLE_DEPRECATED
249  virtual void time_evolving (unsigned int var)
250  {
251  libmesh_deprecated();
252  this->time_evolving(var,1);
253  }
254 #endif
255 
268  virtual void time_evolving (unsigned int var, unsigned int order);
269 
277  bool is_time_evolving (unsigned int var) const
278  {
279  libmesh_assert_less(var,_time_evolving.size());
280  libmesh_assert( _time_evolving[var] == 0 ||
281  _time_evolving[var] == 1 ||
282  _time_evolving[var] == 2 );
283  return _time_evolving[var];
284  }
285 
294  virtual bool eulerian_residual (bool request_jacobian,
295  DiffContext &) {
296  return request_jacobian;
297  }
298 
318  virtual bool mass_residual (bool request_jacobian,
319  DiffContext &) {
320  return request_jacobian;
321  }
322 
335  virtual bool side_mass_residual (bool request_jacobian,
336  DiffContext &) {
337  return request_jacobian;
338  }
339 
356  virtual bool nonlocal_mass_residual (bool request_jacobian,
357  DiffContext & c);
358 
374  virtual bool damping_residual (bool request_jacobian,
375  DiffContext &) {
376  return request_jacobian;
377  }
378 
391  virtual bool side_damping_residual (bool request_jacobian,
392  DiffContext &) {
393  return request_jacobian;
394  }
395 
406  virtual bool nonlocal_damping_residual (bool request_jacobian,
407  DiffContext &) {
408  return request_jacobian;
409  }
410 
411  /*
412  * Prepares the result of a build_context() call for use.
413  *
414  * Most FEMSystem-based problems will need to reimplement this in order to
415  * call FE::get_*() as their particular physics requires.
416  */
417  virtual void init_context(DiffContext &) {}
418 
442  virtual void set_mesh_system(System * sys);
443 
449  const System * get_mesh_system() const;
450 
456 
471  virtual void set_mesh_x_var(unsigned int var);
472 
477  unsigned int get_mesh_x_var() const;
478 
483  virtual void set_mesh_y_var(unsigned int var);
484 
489  unsigned int get_mesh_y_var() const;
490 
495  virtual void set_mesh_z_var(unsigned int var);
496 
501  unsigned int get_mesh_z_var() const;
502 
508  bool _eulerian_time_deriv (bool request_jacobian,
509  DiffContext &);
510 
512  { return !_first_order_vars.empty(); }
513 
517  const std::set<unsigned int> & get_first_order_vars() const
518  { return _first_order_vars; }
519 
520  bool is_first_order_var( unsigned int var ) const
521  { return _first_order_vars.find(var) != _first_order_vars.end(); }
522 
523 
525  { return !_second_order_vars.empty(); }
526 
530  const std::set<unsigned int> & get_second_order_vars() const
531  { return _second_order_vars; }
532 
533  bool is_second_order_var( unsigned int var ) const
534  { return _second_order_vars.find(var) != _second_order_vars.end(); }
535 
536 
537 protected:
538 
543 
548 
554  std::vector<unsigned int> _time_evolving;
555 
559  std::set<unsigned int> _first_order_vars;
560 
564  std::set<unsigned int> _second_order_vars;
565 
571  std::map<unsigned int,unsigned int> _second_order_dot_vars;
572 
573 };
574 
575 // ------------------------------------------------------------
576 // DifferentiablePhysics inline methods
577 
578 
579 inline
581 {
582  // For now we assume that we're doing fully coupled mesh motion
583  // if (sys && sys != this)
584  // libmesh_not_implemented();
585 
586  // For the foreseeable future we'll assume that we keep these
587  // Systems in the same EquationSystems
588  // libmesh_assert_equal_to (&this->get_equation_systems(),
589  // &sys->get_equation_systems());
590 
591  // And for the immediate future this code may not even work
592  libmesh_experimental();
593 
594  _mesh_sys = sys;
595 }
596 
597 
598 
599 inline
601 {
602  _mesh_x_var = var;
603 }
604 
605 
606 
607 inline
609 {
610  _mesh_y_var = var;
611 }
612 
613 
614 
615 inline
617 {
618  _mesh_z_var = var;
619 }
620 
621 
622 
623 inline
625 {
626  return _mesh_sys;
627 }
628 
629 inline
631 {
632  return _mesh_sys;
633 }
634 
635 inline
637 {
638  return _mesh_x_var;
639 }
640 
641 inline
643 {
644  return _mesh_y_var;
645 }
646 
647 inline
649 {
650  return _mesh_z_var;
651 }
652 
653 
654 
655 } // namespace libMesh
656 
657 
658 #endif // LIBMESH_DIFF_PHYSICS_H
bool is_second_order_var(unsigned int var) const
Definition: diff_physics.h:533
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:191
unsigned int get_mesh_x_var() const
Definition: diff_physics.h:636
bool is_first_order_var(unsigned int var) const
Definition: diff_physics.h:520
const unsigned int invalid_uint
Definition: libmesh.h:245
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:124
virtual bool eulerian_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:294
virtual bool side_damping_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:391
std::set< unsigned int > _second_order_vars
Definition: diff_physics.h:564
const System * get_mesh_system() const
Definition: diff_physics.h:624
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:209
virtual bool damping_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:374
std::map< unsigned int, unsigned int > _second_order_dot_vars
Definition: diff_physics.h:571
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:227
bool _eulerian_time_deriv(bool request_jacobian, DiffContext &)
virtual void set_mesh_x_var(unsigned int var)
Definition: diff_physics.h:600
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:318
unsigned int get_mesh_y_var() const
Definition: diff_physics.h:642
const std::set< unsigned int > & get_second_order_vars() const
Definition: diff_physics.h:530
const std::set< unsigned int > & get_first_order_vars() const
Definition: diff_physics.h:517
virtual void set_mesh_y_var(unsigned int var)
Definition: diff_physics.h:608
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
virtual void time_evolving(unsigned int var)
Definition: diff_physics.h:249
virtual void init_context(DiffContext &)
Definition: diff_physics.h:417
virtual bool nonlocal_damping_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:406
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
virtual void set_mesh_z_var(unsigned int var)
Definition: diff_physics.h:616
unsigned int get_mesh_z_var() const
Definition: diff_physics.h:648
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
virtual void init_physics(const System &sys)
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:277
virtual bool side_mass_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:335
std::set< unsigned int > _first_order_vars
Definition: diff_physics.h:559
virtual bool side_time_derivative(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:171
virtual void set_mesh_system(System *sys)
Definition: diff_physics.h:580
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:142
std::vector< unsigned int > _time_evolving
Definition: diff_physics.h:554