libMesh
diff_physics.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 
26 // C++ includes
27 #include <map>
28 #include <memory>
29 #include <set>
30 #include <vector>
31 
32 namespace libMesh
33 {
34 
35 // Forward Declarations
36 class System;
37 class DiffContext;
38 
77 {
78 public:
79 
85  compute_internal_sides (false),
86  _mesh_sys (nullptr),
90  {}
91 
95  virtual ~DifferentiablePhysics ();
96 
100  virtual std::unique_ptr<DifferentiablePhysics> clone_physics() = 0;
101 
105  virtual void clear_physics ();
106 
110  virtual void init_physics (const System & sys);
111 
125  virtual bool element_time_derivative (bool request_jacobian,
126  DiffContext &)
127  {
128  return request_jacobian;
129  }
130 
144  virtual bool element_constraint (bool request_jacobian,
145  DiffContext &)
146  {
147  return request_jacobian;
148  }
149 
157 
174  virtual bool side_time_derivative (bool request_jacobian,
175  DiffContext &)
176  {
177  return request_jacobian;
178  }
179 
195  virtual bool side_constraint (bool request_jacobian,
196  DiffContext &)
197  {
198  return request_jacobian;
199  }
200 
214  virtual bool nonlocal_time_derivative (bool request_jacobian,
215  DiffContext &)
216  {
217  return request_jacobian;
218  }
219 
233  virtual bool nonlocal_constraint (bool request_jacobian,
234  DiffContext &)
235  {
236  return request_jacobian;
237  }
238 
251  virtual void time_evolving (unsigned int var, unsigned int order);
252 
260  bool is_time_evolving (unsigned int var) const
261  {
262  libmesh_assert_less(var,_time_evolving.size());
263  libmesh_assert( _time_evolving[var] == 0 ||
264  _time_evolving[var] == 1 ||
265  _time_evolving[var] == 2 );
266  return _time_evolving[var];
267  }
268 
277  virtual bool eulerian_residual (bool request_jacobian,
278  DiffContext &)
279  {
280  return request_jacobian;
281  }
282 
302  virtual bool mass_residual (bool request_jacobian,
303  DiffContext &)
304  {
305  return request_jacobian;
306  }
307 
320  virtual bool side_mass_residual (bool request_jacobian,
321  DiffContext &)
322  {
323  return request_jacobian;
324  }
325 
342  virtual bool nonlocal_mass_residual (bool request_jacobian,
343  DiffContext & c);
344 
360  virtual bool damping_residual (bool request_jacobian,
361  DiffContext &)
362  {
363  return request_jacobian;
364  }
365 
378  virtual bool side_damping_residual (bool request_jacobian,
379  DiffContext &)
380  {
381  return request_jacobian;
382  }
383 
394  virtual bool nonlocal_damping_residual (bool request_jacobian,
395  DiffContext &)
396  {
397  return request_jacobian;
398  }
399 
400  /*
401  * Prepares the result of a build_context() call for use.
402  *
403  * Most FEMSystem-based problems will need to reimplement this in order to
404  * call FE::get_*() as their particular physics requires.
405  */
406  virtual void init_context(DiffContext &) {}
407 
431  virtual void set_mesh_system(System * sys);
432 
438  const System * get_mesh_system() const;
439 
445 
460  virtual void set_mesh_x_var(unsigned int var);
461 
466  unsigned int get_mesh_x_var() const;
467 
472  virtual void set_mesh_y_var(unsigned int var);
473 
478  unsigned int get_mesh_y_var() const;
479 
484  virtual void set_mesh_z_var(unsigned int var);
485 
490  unsigned int get_mesh_z_var() const;
491 
497  bool _eulerian_time_deriv (bool request_jacobian,
498  DiffContext &);
499 
501  { return !_first_order_vars.empty(); }
502 
506  const std::set<unsigned int> & get_first_order_vars() const
507  { return _first_order_vars; }
508 
509  bool is_first_order_var( unsigned int var ) const
510  { return _first_order_vars.find(var) != _first_order_vars.end(); }
511 
512 
514  { return !_second_order_vars.empty(); }
515 
519  const std::set<unsigned int> & get_second_order_vars() const
520  { return _second_order_vars; }
521 
522  bool is_second_order_var( unsigned int var ) const
523  { return _second_order_vars.find(var) != _second_order_vars.end(); }
524 
525 
526 protected:
527 
532 
537 
543  std::vector<unsigned int> _time_evolving;
544 
548  std::set<unsigned int> _first_order_vars;
549 
553  std::set<unsigned int> _second_order_vars;
554 
560  std::map<unsigned int,unsigned int> _second_order_dot_vars;
561 
562 };
563 
564 // ------------------------------------------------------------
565 // DifferentiablePhysics inline methods
566 
567 
568 inline
570 {
571  // For now we assume that we're doing fully coupled mesh motion
572  // if (sys && sys != this)
573  // libmesh_not_implemented();
574 
575  // For the foreseeable future we'll assume that we keep these
576  // Systems in the same EquationSystems
577  // libmesh_assert_equal_to (&this->get_equation_systems(),
578  // &sys->get_equation_systems());
579 
580  // And for the immediate future this code may not even work
581  libmesh_experimental();
582 
583  _mesh_sys = sys;
584 }
585 
586 
587 
588 inline
590 {
591  _mesh_x_var = var;
592 }
593 
594 
595 
596 inline
598 {
599  _mesh_y_var = var;
600 }
601 
602 
603 
604 inline
606 {
607  _mesh_z_var = var;
608 }
609 
610 
611 
612 inline
614 {
615  return _mesh_sys;
616 }
617 
618 inline
620 {
621  return _mesh_sys;
622 }
623 
624 inline
626 {
627  return _mesh_x_var;
628 }
629 
630 inline
632 {
633  return _mesh_y_var;
634 }
635 
636 inline
638 {
639  return _mesh_z_var;
640 }
641 
642 
643 
644 } // namespace libMesh
645 
646 
647 #endif // LIBMESH_DIFF_PHYSICS_H
bool is_second_order_var(unsigned int var) const
Definition: diff_physics.h:522
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on side of elem to elem_residual.
Definition: diff_physics.h:195
unsigned int get_mesh_x_var() const
Definition: diff_physics.h:625
bool is_first_order_var(unsigned int var) const
Definition: diff_physics.h:509
virtual void clear_physics()
Clear any data structures associated with the physics.
Definition: diff_physics.C:28
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:536
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:286
virtual ~DifferentiablePhysics()
Destructor.
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on elem to elem_residual.
Definition: diff_physics.h:125
virtual bool eulerian_residual(bool request_jacobian, DiffContext &)
Adds a pseudo-convection contribution on elem to elem_residual, if the nodes of elem are being transl...
Definition: diff_physics.h:277
virtual bool side_damping_residual(bool request_jacobian, DiffContext &)
Subtracts a damping vector contribution on side of elem from elem_residual.
Definition: diff_physics.h:378
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:553
const System * get_mesh_system() const
Definition: diff_physics.h:613
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Adds any nonlocal time derivative contributions (e.g.
Definition: diff_physics.h:214
virtual bool damping_residual(bool request_jacobian, DiffContext &)
Subtracts a damping vector contribution on elem from elem_residual.
Definition: diff_physics.h:360
The libMesh namespace provides an interface to certain functionality in the library.
virtual void time_evolving(unsigned int var, unsigned int order)
Tells the DiffSystem that variable var is evolving with respect to time.
Definition: diff_physics.C:41
std::map< unsigned int, unsigned int > _second_order_dot_vars
If the user adds any second order variables, then we need to also cache the map to their correspondin...
Definition: diff_physics.h:560
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Adds any nonlocal constraint contributions (e.g.
Definition: diff_physics.h:233
bool _eulerian_time_deriv(bool request_jacobian, DiffContext &)
This method simply combines element_time_derivative() and eulerian_residual(), which makes its addres...
Definition: diff_physics.C:96
virtual void set_mesh_x_var(unsigned int var)
Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the x...
Definition: diff_physics.h:589
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on elem from elem_residual.
Definition: diff_physics.h:302
virtual void init_physics(const System &sys)
Initialize any data structures associated with the physics.
Definition: diff_physics.C:35
unsigned int get_mesh_y_var() const
Definition: diff_physics.h:631
const std::set< unsigned int > & get_second_order_vars() const
Definition: diff_physics.h:519
const std::set< unsigned int > & get_first_order_vars() const
Definition: diff_physics.h:506
virtual void set_mesh_y_var(unsigned int var)
Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the y...
Definition: diff_physics.h:597
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:531
libmesh_assert(ctx)
DifferentiablePhysics()
Constructor.
Definition: diff_physics.h:84
This class provides a specific system class.
Definition: diff_physics.h:76
virtual void init_context(DiffContext &)
Definition: diff_physics.h:406
bool compute_internal_sides
compute_internal_sides is false by default, indicating that side_* computations will only be done on ...
Definition: diff_physics.h:156
virtual bool nonlocal_damping_residual(bool request_jacobian, DiffContext &)
Subtracts any nonlocal damping vector contributions (e.g.
Definition: diff_physics.h:394
virtual void set_mesh_z_var(unsigned int var)
Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the z...
Definition: diff_physics.h:605
unsigned int get_mesh_z_var() const
Definition: diff_physics.h:637
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
Copy of this object.
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:260
virtual bool side_mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on side of elem from elem_residual.
Definition: diff_physics.h:320
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:548
virtual bool side_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on side of elem to elem_residual.
Definition: diff_physics.h:174
virtual void set_mesh_system(System *sys)
Tells the DifferentiablePhysics that system sys contains the isoparametric Lagrangian variables which...
Definition: diff_physics.h:569
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on elem to elem_residual.
Definition: diff_physics.h:144
std::vector< unsigned int > _time_evolving
Stores unsigned int to tell us which variables are evolving as first order in time (1)...
Definition: diff_physics.h:543
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
Subtracts any nonlocal mass vector contributions (e.g.
Definition: diff_physics.C:57