$extrastylesheet
#include <diff_physics.h>

Public Member Functions | |
| DifferentiablePhysics () | |
| virtual | ~DifferentiablePhysics () |
| virtual UniquePtr < DifferentiablePhysics > | clone_physics ()=0 |
| virtual void | clear_physics () |
| virtual void | init_physics (const System &sys) |
| virtual bool | element_time_derivative (bool request_jacobian, DiffContext &) |
| virtual bool | element_constraint (bool request_jacobian, DiffContext &) |
| virtual bool | side_time_derivative (bool request_jacobian, DiffContext &) |
| virtual bool | side_constraint (bool request_jacobian, DiffContext &) |
| virtual bool | nonlocal_time_derivative (bool request_jacobian, DiffContext &) |
| virtual bool | nonlocal_constraint (bool request_jacobian, DiffContext &) |
| virtual void | time_evolving (unsigned int var) |
| bool | is_time_evolving (unsigned int var) const |
| virtual bool | eulerian_residual (bool request_jacobian, DiffContext &) |
| virtual bool | mass_residual (bool request_jacobian, DiffContext &) |
| virtual bool | side_mass_residual (bool request_jacobian, DiffContext &) |
| virtual bool | nonlocal_mass_residual (bool request_jacobian, DiffContext &c) |
| virtual void | init_context (DiffContext &) |
| virtual void | set_mesh_system (System *sys) |
| const System * | get_mesh_system () const |
| System * | get_mesh_system () |
| virtual void | set_mesh_x_var (unsigned int var) |
| unsigned int | get_mesh_x_var () const |
| virtual void | set_mesh_y_var (unsigned int var) |
| unsigned int | get_mesh_y_var () const |
| virtual void | set_mesh_z_var (unsigned int var) |
| unsigned int | get_mesh_z_var () const |
| bool | _eulerian_time_deriv (bool request_jacobian, DiffContext &) |
Public Attributes | |
| bool | compute_internal_sides |
Protected Attributes | |
| System * | _mesh_sys |
| unsigned int | _mesh_x_var |
| unsigned int | _mesh_y_var |
| unsigned int | _mesh_z_var |
| std::vector< bool > | _time_evolving |
This class provides a specific system class. It aims to generalize any system, linear or nonlinear, which provides both a residual and a Jacobian.
Definition at line 52 of file diff_physics.h.
Constructor. Optionally initializes required data structures.
Definition at line 60 of file diff_physics.h.
:
compute_internal_sides (false),
_mesh_sys (NULL),
_mesh_x_var (libMesh::invalid_uint),
_mesh_y_var (libMesh::invalid_uint),
_mesh_z_var (libMesh::invalid_uint)
{}
| virtual libMesh::DifferentiablePhysics::~DifferentiablePhysics | ( | ) | [virtual] |
Destructor.
| bool libMesh::DifferentiablePhysics::_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.
Referenced by libMesh::EulerSolver::element_residual(), and libMesh::Euler2Solver::element_residual().
| virtual void libMesh::DifferentiablePhysics::clear_physics | ( | ) | [virtual] |
Clear any data structures associated with the physics.
Referenced by libMesh::DifferentiableSystem::clear().
| 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().
| virtual bool libMesh::DifferentiablePhysics::element_constraint | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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().
Definition at line 119 of file diff_physics.h.
Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::SteadySolver::element_residual(), and libMesh::EigenTimeSolver::element_residual().
{
return request_jacobian;
}
| virtual bool libMesh::DifferentiablePhysics::element_time_derivative | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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().
Definition at line 101 of file diff_physics.h.
Referenced by libMesh::SteadySolver::element_residual(), and libMesh::EigenTimeSolver::element_residual().
{
return request_jacobian;
}
| virtual bool libMesh::DifferentiablePhysics::eulerian_residual | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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.
Definition at line 246 of file diff_physics.h.
{
return request_jacobian;
}
| const System * libMesh::DifferentiablePhysics::get_mesh_system | ( | ) | const [inline] |
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 477 of file diff_physics.h.
References _mesh_sys.
Referenced by libMesh::FEMSystem::build_context().
{
return _mesh_sys;
}
| System * libMesh::DifferentiablePhysics::get_mesh_system | ( | ) | [inline] |
Returns a reference to the system with variables corresponding to mesh nodal coordinates, or NULL if the mesh is fixed.
Definition at line 483 of file diff_physics.h.
References _mesh_sys.
{
return _mesh_sys;
}
| unsigned int libMesh::DifferentiablePhysics::get_mesh_x_var | ( | ) | const [inline] |
Returns the variable number corresponding to the mesh x coordinate. Useful for ALE calculations.
Definition at line 489 of file diff_physics.h.
References _mesh_x_var.
Referenced by libMesh::FEMSystem::build_context().
{
return _mesh_x_var;
}
| unsigned int libMesh::DifferentiablePhysics::get_mesh_y_var | ( | ) | const [inline] |
Returns the variable number corresponding to the mesh y coordinate. Useful for ALE calculations.
Definition at line 495 of file diff_physics.h.
References _mesh_y_var.
Referenced by libMesh::FEMSystem::build_context().
{
return _mesh_y_var;
}
| unsigned int libMesh::DifferentiablePhysics::get_mesh_z_var | ( | ) | const [inline] |
Returns the variable number corresponding to the mesh z coordinate. Useful for ALE calculations.
Definition at line 501 of file diff_physics.h.
References _mesh_z_var.
Referenced by libMesh::FEMSystem::build_context().
{
return _mesh_z_var;
}
| virtual void libMesh::DifferentiablePhysics::init_context | ( | DiffContext & | ) | [inline, virtual] |
| virtual void libMesh::DifferentiablePhysics::init_physics | ( | const System & | sys | ) | [virtual] |
Initialize any data structures associated with the physics.
Referenced by libMesh::DifferentiableSystem::attach_physics(), and libMesh::DifferentiableSystem::init_data().
| bool libMesh::DifferentiablePhysics::is_time_evolving | ( | unsigned int | var | ) | const [inline] |
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 234 of file diff_physics.h.
References _time_evolving.
Referenced by libMesh::FEMSystem::init_context().
{
return _time_evolving[var];
}
| virtual bool libMesh::DifferentiablePhysics::mass_residual | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, virtual] |
Subtracts a mass vector contribution on 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.
Many 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 will need to reimplement this themselves.
Reimplemented in libMesh::FEMPhysics.
Definition at line 265 of file diff_physics.h.
Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), and libMesh::EigenTimeSolver::element_residual().
{
return request_jacobian;
}
| virtual bool libMesh::DifferentiablePhysics::nonlocal_constraint | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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-tranient equations have been added.
Definition at line 204 of file diff_physics.h.
Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::SteadySolver::nonlocal_residual(), and libMesh::EigenTimeSolver::nonlocal_residual().
{
return request_jacobian;
}
| virtual 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.
Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), and libMesh::EigenTimeSolver::nonlocal_residual().
| virtual bool libMesh::DifferentiablePhysics::nonlocal_time_derivative | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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 186 of file diff_physics.h.
Referenced by libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::SteadySolver::nonlocal_residual(), and libMesh::EigenTimeSolver::nonlocal_residual().
{
return request_jacobian;
}
| void libMesh::DifferentiablePhysics::set_mesh_system | ( | System * | sys | ) | [inline, 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 433 of file diff_physics.h.
References _mesh_sys, and libMesh::sys.
{
// For now we assume that we're doing fully coupled mesh motion
// if (sys && sys != this)
// libmesh_not_implemented();
// For the foreseeable future we'll assume that we keep these
// Systems in the same EquationSystems
// libmesh_assert_equal_to (&this->get_equation_systems(),
// &sys->get_equation_systems());
// And for the immediate future this code may not even work
libmesh_experimental();
_mesh_sys = sys;
}
| void libMesh::DifferentiablePhysics::set_mesh_x_var | ( | unsigned int | var | ) | [inline, 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 453 of file diff_physics.h.
References _mesh_x_var.
{
_mesh_x_var = var;
}
| void libMesh::DifferentiablePhysics::set_mesh_y_var | ( | unsigned int | var | ) | [inline, 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 461 of file diff_physics.h.
References _mesh_y_var.
{
_mesh_y_var = var;
}
| void libMesh::DifferentiablePhysics::set_mesh_z_var | ( | unsigned int | var | ) | [inline, 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 469 of file diff_physics.h.
References _mesh_z_var.
{
_mesh_z_var = var;
}
| virtual bool libMesh::DifferentiablePhysics::side_constraint | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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().
Definition at line 168 of file diff_physics.h.
Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::SteadySolver::side_residual(), and libMesh::EigenTimeSolver::side_residual().
{
return request_jacobian;
}
| virtual bool libMesh::DifferentiablePhysics::side_mass_residual | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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 282 of file diff_physics.h.
Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), and libMesh::EigenTimeSolver::side_residual().
{
return request_jacobian;
}
| virtual bool libMesh::DifferentiablePhysics::side_time_derivative | ( | bool | request_jacobian, |
| DiffContext & | |||
| ) | [inline, 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().
Definition at line 148 of file diff_physics.h.
Referenced by libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::SteadySolver::side_residual(), and libMesh::EigenTimeSolver::side_residual().
{
return request_jacobian;
}
| virtual void libMesh::DifferentiablePhysics::time_evolving | ( | unsigned int | var | ) | [inline, 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 reimplment this function; however any system which reimplements mass_residual() may have to reimplement time_evolving() to prepare data structures.
Definition at line 221 of file diff_physics.h.
References _time_evolving.
{
if (_time_evolving.size() <= var)
_time_evolving.resize(var+1, false);
_time_evolving[var] = true;
}
System* libMesh::DifferentiablePhysics::_mesh_sys [protected] |
System from which to acquire moving mesh information
Definition at line 412 of file diff_physics.h.
Referenced by get_mesh_system(), libMesh::FEMSystem::mesh_position_get(), libMesh::FEMSystem::mesh_position_set(), libMesh::FEMSystem::numerical_jacobian(), and set_mesh_system().
unsigned int libMesh::DifferentiablePhysics::_mesh_x_var [protected] |
Variables from which to acquire moving mesh information
Definition at line 417 of file diff_physics.h.
Referenced by 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] |
Definition at line 417 of file diff_physics.h.
Referenced by get_mesh_y_var(), libMesh::FEMSystem::mesh_position_get(), libMesh::FEMSystem::numerical_jacobian(), and set_mesh_y_var().
unsigned int libMesh::DifferentiablePhysics::_mesh_z_var [protected] |
Definition at line 417 of file diff_physics.h.
Referenced by get_mesh_z_var(), libMesh::FEMSystem::mesh_position_get(), libMesh::FEMSystem::numerical_jacobian(), and set_mesh_z_var().
std::vector<bool> libMesh::DifferentiablePhysics::_time_evolving [protected] |
Stores bools to tell us which variables are evolving in time and which are just constraints
Definition at line 423 of file diff_physics.h.
Referenced by is_time_evolving(), and time_evolving().
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 130 of file diff_physics.h.