libMesh
Public Member Functions | Public Attributes | Protected Attributes | List of all members
libMesh::DifferentiablePhysics Class Referenceabstract

This class provides a specific system class. More...

#include <diff_physics.h>

Inheritance diagram for libMesh::DifferentiablePhysics:
[legend]

Public Member Functions

 DifferentiablePhysics ()
 Constructor. More...
 
virtual ~DifferentiablePhysics ()
 Destructor. More...
 
virtual UniquePtr< DifferentiablePhysicsclone_physics ()=0
 Copy of this object. More...
 
virtual void clear_physics ()
 Clear any data structures associated with the physics. More...
 
virtual void init_physics (const System &sys)
 Initialize any data structures associated with the physics. More...
 
virtual bool element_time_derivative (bool request_jacobian, DiffContext &)
 Adds the time derivative contribution on elem to elem_residual. More...
 
virtual bool element_constraint (bool request_jacobian, DiffContext &)
 Adds the constraint contribution on elem to elem_residual. More...
 
virtual bool side_time_derivative (bool request_jacobian, DiffContext &)
 Adds the time derivative contribution on side of elem to elem_residual. More...
 
virtual bool side_constraint (bool request_jacobian, DiffContext &)
 Adds the constraint contribution on side of elem to elem_residual. More...
 
virtual bool nonlocal_time_derivative (bool request_jacobian, DiffContext &)
 Adds any nonlocal time derivative contributions (e.g. More...
 
virtual bool nonlocal_constraint (bool request_jacobian, DiffContext &)
 Adds any nonlocal constraint contributions (e.g. More...
 
virtual void time_evolving (unsigned int var)
 Tells the DiffSystem that variable var is evolving with respect to time. More...
 
virtual void time_evolving (unsigned int var, unsigned int order)
 Tells the DiffSystem that variable var is evolving with respect to time. More...
 
bool is_time_evolving (unsigned int var) const
 
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 translated by a moving mesh. More...
 
virtual bool mass_residual (bool request_jacobian, DiffContext &)
 Subtracts a mass vector contribution on elem from elem_residual. More...
 
virtual bool side_mass_residual (bool request_jacobian, DiffContext &)
 Subtracts a mass vector contribution on side of elem from elem_residual. More...
 
virtual bool nonlocal_mass_residual (bool request_jacobian, DiffContext &c)
 Subtracts any nonlocal mass vector contributions (e.g. More...
 
virtual bool damping_residual (bool request_jacobian, DiffContext &)
 Subtracts a damping vector contribution on elem from elem_residual. More...
 
virtual bool side_damping_residual (bool request_jacobian, DiffContext &)
 Subtracts a damping vector contribution on side of elem from elem_residual. More...
 
virtual bool nonlocal_damping_residual (bool request_jacobian, DiffContext &)
 Subtracts any nonlocal damping vector contributions (e.g. More...
 
virtual void init_context (DiffContext &)
 
virtual void set_mesh_system (System *sys)
 Tells the DifferentiablePhysics that system sys contains the isoparametric Lagrangian variables which correspond to the coordinates of mesh nodes, in problems where the mesh itself is expected to move in time. More...
 
const Systemget_mesh_system () const
 
Systemget_mesh_system ()
 
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 coordinate of mesh nodes, in problems where the mesh itself is expected to move in time. More...
 
unsigned int get_mesh_x_var () const
 
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 coordinate of mesh nodes. More...
 
unsigned int get_mesh_y_var () const
 
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 coordinate of mesh nodes. More...
 
unsigned int get_mesh_z_var () const
 
bool _eulerian_time_deriv (bool request_jacobian, DiffContext &)
 This method simply combines element_time_derivative() and eulerian_residual(), which makes its address useful as a pointer-to-member-function when refactoring. More...
 
bool have_first_order_vars () const
 
const std::set< unsigned int > & get_first_order_vars () const
 
bool is_first_order_var (unsigned int var) const
 
bool have_second_order_vars () const
 
const std::set< unsigned int > & get_second_order_vars () const
 
bool is_second_order_var (unsigned int var) const
 

Public Attributes

bool compute_internal_sides
 compute_internal_sides is false by default, indicating that side_* computations will only be done on boundary sides. More...
 

Protected Attributes

System_mesh_sys
 System from which to acquire moving mesh information. More...
 
unsigned int _mesh_x_var
 Variables from which to acquire moving mesh information. More...
 
unsigned int _mesh_y_var
 
unsigned int _mesh_z_var
 
std::vector< unsigned int_time_evolving
 Stores unsigned int to tell us which variables are evolving as first order in time (1), second order in time (2), or are not time evolving (0). More...
 
std::set< unsigned int_first_order_vars
 Variable indices for those variables that are first order in time. More...
 
std::set< unsigned int_second_order_vars
 Variable indices for those variables that are second order in time. More...
 
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 corresponding dot variable that will be added by this TimeSolver class. More...
 

Detailed Description

This class provides a specific system class.

It aims to generalize any system, linear or nonlinear, which provides both a residual and a Jacobian. For first order (in time) systems, the (nonlinear) residual computed at each time step is

\[ F(u) + G(u) - M(u,\dot{u})\dot{u} = 0 \]

for unsteady TimeSolver and $ F(u) = 0$ for steady TimeSolver. $F(u)$ is computed by element/side_time_derivative, $ G(u) $ is computed using element/side_constraint, and $ M(u,\dot{u})\dot{u} $ is computed using the mass_residual methods.

For second order (in time) systems, the (nonlinear) residual computed at each time step is

\[ M(u,\ddot{u})\ddot{u} + C(u,\dot{u})\dot{u} + F(u) + G(u) = 0 \]

for unsteady TimeSolver and $ F(u) = 0$ for steady TimeSolver. $F(u)$ is computed by element/side_time_derivative, $G(u)$ is computed using element/side_constraint, $C(u,\dot{u})\dot{u}$ is computed using the damping_residual methods and $ -M(u,\ddot{u})\ddot{u}$ is computed using the mass_residual methods. This is the sign convention used by the default implementation; if the method is overridden, the user can choose any self-consistent sign convention they wish.

FEMContext provides methods for querying values of the solution ${u}$, its "rate" $\dot{u}$ and its "acceleration" $\ddot{u}$. Furthermore, derivatives of each of these w.r.t the nonlinear iteration unknown (e.g. in EulerSolver, the solution at the next time step $ u_{n+1} $) are provided through DiffContext::get_elem_solution_derivative(), DiffContext::get_elem_solution_rate_derivative(), and DiffContext::get_elem_solution_accel_derivative(). The should be incorporated into the Jacobian evaluations, if the Jacobian is being provided.

Author
Roy H. Stogner
Date
2006

Definition at line 74 of file diff_physics.h.

Constructor & Destructor Documentation

libMesh::DifferentiablePhysics::DifferentiablePhysics ( )

Constructor.

Optionally initializes required data structures.

Definition at line 82 of file diff_physics.h.

References clear_physics(), clone_physics(), init_physics(), libMesh::sys, and ~DifferentiablePhysics().

82  :
83  compute_internal_sides (false),
88  {}
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:546
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:184
const class libmesh_nullptr_t libmesh_nullptr
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:541
bool compute_internal_sides
compute_internal_sides is false by default, indicating that side_* computations will only be done on ...
Definition: diff_physics.h:152
libMesh::DifferentiablePhysics::~DifferentiablePhysics ( )
virtual

Destructor.

Definition at line 26 of file diff_physics.C.

References clear_physics().

Referenced by DifferentiablePhysics().

27 {
29 }
virtual void clear_physics()
Clear any data structures associated with the physics.
Definition: diff_physics.C:33

Member Function Documentation

bool libMesh::DifferentiablePhysics::_eulerian_time_deriv ( bool  request_jacobian,
DiffContext context 
)

This method simply combines element_time_derivative() and eulerian_residual(), which makes its address useful as a pointer-to-member-function when refactoring.

Definition at line 102 of file diff_physics.C.

References element_time_derivative(), and eulerian_residual().

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::NewmarkSolver::element_residual(), and init_context().

104 {
105  // For any problem we need time derivative terms
106  request_jacobian =
107  this->element_time_derivative(request_jacobian, context);
108 
109  // For a moving mesh problem we may need the pseudoconvection term too
110  return this->eulerian_residual(request_jacobian, context) &&
111  request_jacobian;
112 }
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on elem to elem_residual.
Definition: diff_physics.h:123
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:293
void libMesh::DifferentiablePhysics::clear_physics ( )
virtual

Clear any data structures associated with the physics.

Definition at line 33 of file diff_physics.C.

References _time_evolving.

Referenced by libMesh::DifferentiableSystem::clear(), DifferentiablePhysics(), and ~DifferentiablePhysics().

34 {
35  _time_evolving.resize(0);
36 }
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:553
virtual UniquePtr<DifferentiablePhysics> libMesh::DifferentiablePhysics::clone_physics ( )
pure virtual

Copy of this object.

User should override to copy any needed state.

Implemented in libMesh::DifferentiableSystem.

Referenced by libMesh::DifferentiableSystem::attach_physics(), and DifferentiablePhysics().

virtual bool libMesh::DifferentiablePhysics::damping_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Subtracts a damping vector contribution on elem from elem_residual.

This method is not used in first-order-in-time problems. For second-order-in-time problems, this is the $ C(u,\ddot{u})\ddot{u} $ term. This method is only called for UnsteadySolver-based TimeSolvers.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

If the problem has no damping, the default "do-nothing" is correct. Otherwise, this must be reimplemented.

Reimplemented in SecondOrderScalarSystemFirstOrderTimeSolverBase, and SecondOrderScalarSystemSecondOrderTimeSolverBase.

Definition at line 373 of file diff_physics.h.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), and libMesh::NewmarkSolver::element_residual().

374  {
375  return request_jacobian;
376  }
virtual bool libMesh::DifferentiablePhysics::element_constraint ( bool  request_jacobian,
DiffContext  
)
virtual

Adds the constraint contribution on elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users may need to reimplement this for their particular PDE.

To implement the constraint 0 = G(u), the user should examine u = elem_solution and add (G(u), phi_i) to elem_residual in elem_constraint().

Reimplemented in CoupledSystem, and NavierSystem.

Definition at line 141 of file diff_physics.h.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::SteadySolver::element_residual(), libMesh::EigenTimeSolver::element_residual(), and libMesh::NewmarkSolver::element_residual().

142  {
143  return request_jacobian;
144  }
virtual bool libMesh::DifferentiablePhysics::element_time_derivative ( bool  request_jacobian,
DiffContext  
)
virtual

Adds the time derivative contribution on elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users need to reimplement this for their particular PDE.

To implement the physics model du/dt = F(u), the user should examine u = elem_solution and add (F(u), phi_i) to elem_residual in elem_time_derivative().

Reimplemented in SecondOrderScalarSystemFirstOrderTimeSolverBase, FirstOrderScalarSystemBase, HeatSystem, CoupledSystem, ElasticitySystem, LaplaceSystem, CurlCurlSystem, LaplaceSystem, PoissonSystem, LaplaceSystem, LaplaceSystem, CurlCurlSystem, L2System, SolidSystem, NavierSystem, and HeatSystem.

Definition at line 123 of file diff_physics.h.

Referenced by _eulerian_time_deriv(), libMesh::SteadySolver::element_residual(), and libMesh::EigenTimeSolver::element_residual().

124  {
125  return request_jacobian;
126  }
virtual bool libMesh::DifferentiablePhysics::eulerian_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Adds a pseudo-convection contribution on elem to elem_residual, if the nodes of elem are being translated by a moving mesh.

The library provides a basic implementation in FEMPhysics::eulerian_residual()

Reimplemented in libMesh::FEMPhysics, and SolidSystem.

Definition at line 293 of file diff_physics.h.

Referenced by _eulerian_time_deriv().

294  {
295  return request_jacobian;
296  }
const std::set<unsigned int>& libMesh::DifferentiablePhysics::get_first_order_vars ( ) const
Returns
The set of first order in time variable indices. May be empty.

Definition at line 516 of file diff_physics.h.

References _first_order_vars.

Referenced by libMesh::DifferentiableSystem::have_first_order_scalar_vars().

517  { return _first_order_vars; }
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:558
const System * libMesh::DifferentiablePhysics::get_mesh_system ( ) const
Returns
A const reference to the system with variables corresponding to mesh nodal coordinates, or NULL if the mesh is fixed. Useful for ALE calculations.

Definition at line 623 of file diff_physics.h.

References _mesh_sys.

Referenced by libMesh::FEMSystem::build_context(), and init_context().

624 {
625  return _mesh_sys;
626 }
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:541
System * libMesh::DifferentiablePhysics::get_mesh_system ( )
Returns
A reference to the system with variables corresponding to mesh nodal coordinates, or NULL if the mesh is fixed.

Definition at line 629 of file diff_physics.h.

References _mesh_sys.

630 {
631  return _mesh_sys;
632 }
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:541
unsigned int libMesh::DifferentiablePhysics::get_mesh_x_var ( ) const
Returns
The variable number corresponding to the mesh x coordinate. Useful for ALE calculations.

Definition at line 635 of file diff_physics.h.

References _mesh_x_var.

Referenced by libMesh::FEMSystem::build_context(), and init_context().

636 {
637  return _mesh_x_var;
638 }
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:546
unsigned int libMesh::DifferentiablePhysics::get_mesh_y_var ( ) const
Returns
The variable number corresponding to the mesh y coordinate. Useful for ALE calculations.

Definition at line 641 of file diff_physics.h.

References _mesh_y_var.

Referenced by libMesh::FEMSystem::build_context(), and init_context().

642 {
643  return _mesh_y_var;
644 }
unsigned int libMesh::DifferentiablePhysics::get_mesh_z_var ( ) const
Returns
The variable number corresponding to the mesh z coordinate. Useful for ALE calculations.

Definition at line 647 of file diff_physics.h.

References _mesh_z_var.

Referenced by libMesh::FEMSystem::build_context(), and init_context().

648 {
649  return _mesh_z_var;
650 }
const std::set<unsigned int>& libMesh::DifferentiablePhysics::get_second_order_vars ( ) const
Returns
The set of second order in time variable indices. May be empty.

Definition at line 529 of file diff_physics.h.

References _second_order_vars.

Referenced by libMesh::DifferentiableSystem::add_second_order_dot_vars(), libMesh::DiffContext::DiffContext(), libMesh::Euler2Solver::element_residual(), libMesh::DifferentiableSystem::have_second_order_scalar_vars(), and libMesh::FEMContext::pre_fe_reinit().

530  { return _second_order_vars; }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:563
bool libMesh::DifferentiablePhysics::have_first_order_vars ( ) const

Definition at line 510 of file diff_physics.h.

References _first_order_vars.

Referenced by libMesh::DifferentiableSystem::have_first_order_scalar_vars().

511  { return !_first_order_vars.empty(); }
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:558
bool libMesh::DifferentiablePhysics::have_second_order_vars ( ) const

Definition at line 523 of file diff_physics.h.

References _second_order_vars.

Referenced by libMesh::EulerSolver::element_residual(), and libMesh::DifferentiableSystem::have_second_order_scalar_vars().

524  { return !_second_order_vars.empty(); }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:563
virtual void libMesh::DifferentiablePhysics::init_context ( DiffContext )
virtual
void libMesh::DifferentiablePhysics::init_physics ( const System sys)
virtual

Initialize any data structures associated with the physics.

Definition at line 40 of file diff_physics.C.

References _time_evolving, and libMesh::System::n_vars().

Referenced by libMesh::DifferentiableSystem::attach_physics(), DifferentiablePhysics(), and libMesh::DifferentiableSystem::init_data().

41 {
42  // give us flags for every variable that might be time evolving
43  _time_evolving.resize(sys.n_vars(), false);
44 }
ImplicitSystem & sys
unsigned int n_vars() const
Definition: system.h:2086
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:553
bool libMesh::DifferentiablePhysics::is_first_order_var ( unsigned int  var) const

Definition at line 519 of file diff_physics.h.

References _first_order_vars.

520  { return _first_order_vars.find(var) != _first_order_vars.end(); }
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:558
bool libMesh::DifferentiablePhysics::is_second_order_var ( unsigned int  var) const

Definition at line 532 of file diff_physics.h.

References _second_order_vars.

Referenced by libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns().

533  { return _second_order_vars.find(var) != _second_order_vars.end(); }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:563
bool libMesh::DifferentiablePhysics::is_time_evolving ( unsigned int  var) const
Returns
true iff variable var is evolving with respect to time. In general, the user's init() function should have set time_evolving() for any variables which behave like du/dt = F(u), and should not call time_evolving() for any variables which behave like 0 = G(u).

Definition at line 276 of file diff_physics.h.

References _time_evolving, and libMesh::libmesh_assert().

Referenced by libMesh::FEMPhysics::eulerian_residual(), libMesh::FEMSystem::init_context(), libMesh::FEMPhysics::mass_residual(), and nonlocal_mass_residual().

277  {
278  libmesh_assert_less(var,_time_evolving.size());
279  libmesh_assert( _time_evolving[var] == 0 ||
280  _time_evolving[var] == 1 ||
281  _time_evolving[var] == 2 );
282  return _time_evolving[var];
283  }
libmesh_assert(j)
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:553
virtual bool libMesh::DifferentiablePhysics::mass_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Subtracts a mass vector contribution on elem from elem_residual.

For first-order-in-time problems, this is the $ M(u,\dot{u})\dot{u} $ term. For second-order-in-time problems, this is the $ M(u,\ddot{u})\ddot{u} $ term. This method is only called for UnsteadySolver-based TimeSolvers.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Many first-order-in-time problems can use the reimplementation in FEMPhysics::mass_residual which subtracts (du/dt,v) for each transient variable u; users with more complicated transient problems or second-order-in-time problems will need to reimplement this themselves.

Reimplemented in SecondOrderScalarSystemFirstOrderTimeSolverBase, SecondOrderScalarSystemSecondOrderTimeSolverBase, FirstOrderScalarSystemBase, libMesh::FEMPhysics, ElasticitySystem, and NavierSystem.

Definition at line 317 of file diff_physics.h.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::NewmarkSolver::element_residual(), and libMesh::EigenTimeSolver::element_residual().

318  {
319  return request_jacobian;
320  }
virtual bool libMesh::DifferentiablePhysics::nonlocal_constraint ( bool  request_jacobian,
DiffContext  
)
virtual

Adds any nonlocal constraint contributions (e.g.

some components of constraints in scalar variable equations) to elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Users may need to reimplement this for PDEs on systems to which SCALAR variables with non-transient equations have been added.

Definition at line 226 of file diff_physics.h.

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::SteadySolver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

227  {
228  return request_jacobian;
229  }
virtual bool libMesh::DifferentiablePhysics::nonlocal_damping_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Subtracts any nonlocal damping vector contributions (e.g.

any first time derivative coefficients in scalar variable equations) from elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Definition at line 405 of file diff_physics.h.

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

406  {
407  return request_jacobian;
408  }
bool libMesh::DifferentiablePhysics::nonlocal_mass_residual ( bool  request_jacobian,
DiffContext c 
)
virtual

Subtracts any nonlocal mass vector contributions (e.g.

any time derivative coefficients in scalar variable equations) from elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Many problems can use the reimplementation in FEMPhysics::mass_residual which subtracts (du/dt,v) for each transient scalar variable u; users with more complicated transient scalar variable equations will need to reimplement this themselves.

Definition at line 63 of file diff_physics.C.

References libMesh::DiffContext::elem_solution_rate_derivative, libMesh::FEType::family, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_residual(), libMesh::DiffContext::get_elem_solution(), libMesh::DiffContext::get_system(), is_time_evolving(), libMesh::DiffContext::n_vars(), libMesh::SCALAR, libMesh::Variable::type(), and libMesh::System::variable().

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), libMesh::NewmarkSolver::nonlocal_residual(), and side_mass_residual().

65 {
66  FEMContext & context = cast_ref<FEMContext &>(c);
67 
68  for (unsigned int var = 0; var != context.n_vars(); ++var)
69  {
70  if (!this->is_time_evolving(var))
71  continue;
72 
73  if (c.get_system().variable(var).type().family != SCALAR)
74  continue;
75 
76  const std::vector<dof_id_type> & dof_indices =
77  context.get_dof_indices(var);
78 
79  const unsigned int n_dofs = cast_int<unsigned int>
80  (dof_indices.size());
81 
82  DenseSubVector<Number> & Fs = context.get_elem_residual(var);
83  DenseSubMatrix<Number> & Kss = context.get_elem_jacobian( var, var );
84 
86  context.get_elem_solution(var);
87 
88  for (unsigned int i=0; i != n_dofs; ++i)
89  {
90  Fs(i) -= Us(i);
91 
92  if (request_jacobian)
93  Kss(i,i) -= context.elem_solution_rate_derivative;
94  }
95  }
96 
97  return request_jacobian;
98 }
bool is_time_evolving(unsigned int var) const
Definition: diff_physics.h:276
Defines a dense subvector for use in finite element computations.
virtual bool libMesh::DifferentiablePhysics::nonlocal_time_derivative ( bool  request_jacobian,
DiffContext  
)
virtual

Adds any nonlocal time derivative contributions (e.g.

some components of time derivatives in scalar variable equations) to elem_residual

If this method receives request_jacobian = true, then it should also modify elem_jacobian and return true if possible. If the Jacobian changes have not been computed then the method should return false.

Users may need to reimplement this for PDEs on systems to which SCALAR variables have been added.

Definition at line 208 of file diff_physics.h.

Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::SteadySolver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), and libMesh::NewmarkSolver::nonlocal_residual().

209  {
210  return request_jacobian;
211  }
void libMesh::DifferentiablePhysics::set_mesh_system ( System sys)
virtual

Tells the DifferentiablePhysics that system sys contains the isoparametric Lagrangian variables which correspond to the coordinates of mesh nodes, in problems where the mesh itself is expected to move in time.

The system with mesh coordinate data (which may be this system itself, for fully coupled moving mesh problems) is currently assumed to have new (end of time step) mesh coordinates stored in solution, old (beginning of time step) mesh coordinates stored in _old_nonlinear_solution, and constant velocity motion during each time step.

Activating this function ensures that local (but not neighbor!) element geometry is correctly repositioned when evaluating element residuals.

Currently sys must be *this for a tightly coupled moving mesh problem or NULL to stop mesh movement; loosely coupled moving mesh problems are not implemented.

This code is experimental. "Trust but verify, and not in that order"

Definition at line 579 of file diff_physics.h.

References _mesh_sys, and libMesh::sys.

Referenced by init_context(), and SolidSystem::init_data().

580 {
581  // For now we assume that we're doing fully coupled mesh motion
582  // if (sys && sys != this)
583  // libmesh_not_implemented();
584 
585  // For the foreseeable future we'll assume that we keep these
586  // Systems in the same EquationSystems
587  // libmesh_assert_equal_to (&this->get_equation_systems(),
588  // &sys->get_equation_systems());
589 
590  // And for the immediate future this code may not even work
591  libmesh_experimental();
592 
593  _mesh_sys = sys;
594 }
ImplicitSystem & sys
System * _mesh_sys
System from which to acquire moving mesh information.
Definition: diff_physics.h:541
void libMesh::DifferentiablePhysics::set_mesh_x_var ( unsigned int  var)
virtual

Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the x coordinate of mesh nodes, in problems where the mesh itself is expected to move in time.

The system with mesh coordinate data (which may be this system itself, for fully coupled moving mesh problems) is currently assumed to have new (end of time step) mesh coordinates stored in solution, old (beginning of time step) mesh coordinates stored in _old_nonlinear_solution, and constant velocity motion during each time step.

Activating this function ensures that local (but not neighbor!) element geometry is correctly repositioned when evaluating element residuals.

Definition at line 599 of file diff_physics.h.

References _mesh_x_var.

Referenced by init_context(), and SolidSystem::init_data().

600 {
601  _mesh_x_var = var;
602 }
unsigned int _mesh_x_var
Variables from which to acquire moving mesh information.
Definition: diff_physics.h:546
void libMesh::DifferentiablePhysics::set_mesh_y_var ( unsigned int  var)
virtual

Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the y coordinate of mesh nodes.

Definition at line 607 of file diff_physics.h.

References _mesh_y_var.

Referenced by init_context(), and SolidSystem::init_data().

608 {
609  _mesh_y_var = var;
610 }
void libMesh::DifferentiablePhysics::set_mesh_z_var ( unsigned int  var)
virtual

Tells the DifferentiablePhysics that variable var from the mesh system should be used to update the z coordinate of mesh nodes.

Definition at line 615 of file diff_physics.h.

References _mesh_z_var.

Referenced by init_context(), and SolidSystem::init_data().

616 {
617  _mesh_z_var = var;
618 }
virtual bool libMesh::DifferentiablePhysics::side_constraint ( bool  request_jacobian,
DiffContext  
)
virtual

Adds the constraint contribution on side of elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users may need to reimplement this for their particular PDE depending on the boundary conditions.

To implement a weak form of the constraint 0 = G(u), the user should examine u = elem_solution and add (G(u), phi_i) boundary integral contributions to elem_residual in side_constraint().

Reimplemented in LaplaceSystem, LaplaceSystem, LaplaceSystem, and NavierSystem.

Definition at line 190 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::SteadySolver::side_residual(), libMesh::EigenTimeSolver::side_residual(), and libMesh::NewmarkSolver::side_residual().

191  {
192  return request_jacobian;
193  }
virtual bool libMesh::DifferentiablePhysics::side_damping_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Subtracts a damping vector contribution on side of elem from elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

For most problems, the default implementation of "do nothing" is correct; users with boundary conditions including first time derivatives may need to reimplement this themselves.

Definition at line 390 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), and libMesh::NewmarkSolver::side_residual().

391  {
392  return request_jacobian;
393  }
virtual bool libMesh::DifferentiablePhysics::side_mass_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Subtracts a mass vector contribution on side of elem from elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

For most problems, the default implementation of "do nothing" is correct; users with boundary conditions including time derivatives may need to reimplement this themselves.

Definition at line 334 of file diff_physics.h.

References nonlocal_mass_residual().

Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::EigenTimeSolver::side_residual(), and libMesh::NewmarkSolver::side_residual().

335  {
336  return request_jacobian;
337  }
virtual bool libMesh::DifferentiablePhysics::side_time_derivative ( bool  request_jacobian,
DiffContext  
)
virtual

Adds the time derivative contribution on side of elem to elem_residual.

If this method receives request_jacobian = true, then it should compute elem_jacobian and return true if possible. If elem_jacobian has not been computed then the method should return false.

Users may need to reimplement this for their particular PDE depending on the boundary conditions.

To implement a weak form of the source term du/dt = F(u) on sides, such as might arise in a flux boundary condition, the user should examine u = elem_solution and add (F(u), phi_i) boundary integral contributions to elem_residual in side_constraint().

Reimplemented in ElasticitySystem, CurlCurlSystem, CurlCurlSystem, and SolidSystem.

Definition at line 170 of file diff_physics.h.

Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::SteadySolver::side_residual(), libMesh::EigenTimeSolver::side_residual(), and libMesh::NewmarkSolver::side_residual().

171  {
172  return request_jacobian;
173  }
virtual void libMesh::DifferentiablePhysics::time_evolving ( unsigned int  var)
virtual

Tells the DiffSystem that variable var is evolving with respect to time.

In general, the user's init() function should call time_evolving() for any variables which behave like du/dt = F(u), and should not call time_evolving() for any variables which behave like 0 = G(u).

Most derived systems will not have to reimplement this function; however any system which reimplements mass_residual() may have to reimplement time_evolving() to prepare data structures.

Deprecated:
Instead, use the time_evolving override and specify the order-in-time of the variable, either 1 or 2. This method assumes the variable is first order for backward compatibility.

Definition at line 248 of file diff_physics.h.

Referenced by libMesh::DifferentiableSystem::add_second_order_dot_vars(), CurlCurlSystem::init_data(), and HeatSystem::init_data().

249  {
250  libmesh_deprecated();
251  this->time_evolving(var,1);
252  }
virtual void time_evolving(unsigned int var)
Tells the DiffSystem that variable var is evolving with respect to time.
Definition: diff_physics.h:248
void libMesh::DifferentiablePhysics::time_evolving ( unsigned int  var,
unsigned int  order 
)
virtual

Tells the DiffSystem that variable var is evolving with respect to time.

In general, the user's init() function should call time_evolving() with order 1 for any variables which behave like du/dt = F(u), with order 2 for any variables that behave like d^2u/dt^2 = F(u), and should not call time_evolving() for any variables which behave like 0 = G(u).

Most derived systems will not have to reimplement this function; however any system which reimplements mass_residual() may have to reimplement time_evolving() to prepare data structures.

Definition at line 46 of file diff_physics.C.

References _first_order_vars, _second_order_vars, and _time_evolving.

48 {
49  if (order != 1 && order != 2)
50  libmesh_error_msg("Input order must be 1 or 2!");
51 
52  if (_time_evolving.size() <= var)
53  _time_evolving.resize(var+1, 0);
54 
55  _time_evolving[var] = order;
56 
57  if (order == 1)
58  _first_order_vars.insert(var);
59  else
60  _second_order_vars.insert(var);
61 }
std::set< unsigned int > _second_order_vars
Variable indices for those variables that are second order in time.
Definition: diff_physics.h:563
std::set< unsigned int > _first_order_vars
Variable indices for those variables that are first order in time.
Definition: diff_physics.h:558
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:553

Member Data Documentation

std::set<unsigned int> libMesh::DifferentiablePhysics::_first_order_vars
protected

Variable indices for those variables that are first order in time.

Definition at line 558 of file diff_physics.h.

Referenced by get_first_order_vars(), have_first_order_vars(), is_first_order_var(), and time_evolving().

System* libMesh::DifferentiablePhysics::_mesh_sys
protected
unsigned int libMesh::DifferentiablePhysics::_mesh_x_var
protected

Variables from which to acquire moving mesh information.

Definition at line 546 of file diff_physics.h.

Referenced by libMesh::FEMPhysics::eulerian_residual(), get_mesh_x_var(), libMesh::FEMSystem::mesh_position_get(), libMesh::FEMSystem::numerical_jacobian(), and set_mesh_x_var().

unsigned int libMesh::DifferentiablePhysics::_mesh_y_var
protected
unsigned int libMesh::DifferentiablePhysics::_mesh_z_var
protected
std::map<unsigned int,unsigned int> libMesh::DifferentiablePhysics::_second_order_dot_vars
protected

If the user adds any second order variables, then we need to also cache the map to their corresponding dot variable that will be added by this TimeSolver class.

Definition at line 570 of file diff_physics.h.

Referenced by libMesh::DifferentiableSystem::add_second_order_dot_vars(), and libMesh::DifferentiableSystem::get_second_order_dot_var().

std::set<unsigned int> libMesh::DifferentiablePhysics::_second_order_vars
protected

Variable indices for those variables that are second order in time.

Definition at line 563 of file diff_physics.h.

Referenced by get_second_order_vars(), have_second_order_vars(), is_second_order_var(), and time_evolving().

std::vector<unsigned int> libMesh::DifferentiablePhysics::_time_evolving
protected

Stores unsigned int to tell us which variables are evolving as first order in time (1), second order in time (2), or are not time evolving (0).

Definition at line 553 of file diff_physics.h.

Referenced by clear_physics(), init_physics(), is_time_evolving(), and time_evolving().

bool libMesh::DifferentiablePhysics::compute_internal_sides

compute_internal_sides is false by default, indicating that side_* computations will only be done on boundary sides.

If compute_internal_sides is true, computations will be done on sides between elements as well.

Definition at line 152 of file diff_physics.h.


The documentation for this class was generated from the following files: