libMesh::FEMPhysics Class Referenceabstract

#include <fem_physics.h>

Inheritance diagram for libMesh::FEMPhysics:

Public Member Functions

 FEMPhysics ()
 
virtual ~FEMPhysics ()
 
virtual bool eulerian_residual (bool request_jacobian, DiffContext &context)
 
virtual bool mass_residual (bool request_jacobian, DiffContext &)
 
virtual AutoPtr
< 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 void time_evolving (unsigned int var)
 
bool is_time_evolving (unsigned int var) const
 
virtual bool side_mass_residual (bool request_jacobian, DiffContext &)
 
virtual void init_context (DiffContext &)
 
virtual void set_mesh_system (System *sys)
 
const Systemget_mesh_system () const
 
Systemget_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
 

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
 

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.

This class is part of the new DifferentiableSystem framework, which is still experimental. Users of this framework should beware of bugs and future API changes.

Author
Roy H. Stogner 2012

Definition at line 48 of file fem_physics.h.

Constructor & Destructor Documentation

libMesh::FEMPhysics::FEMPhysics ( )
inline

Constructor.

Definition at line 55 of file fem_physics.h.

55  :
57  {}
virtual libMesh::FEMPhysics::~FEMPhysics ( )
inlinevirtual

Destructor.

Definition at line 62 of file fem_physics.h.

62 {}

Member Function Documentation

virtual void libMesh::DifferentiablePhysics::clear_physics ( )
virtualinherited

Clear any data structures associated with the physics.

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

virtual AutoPtr<DifferentiablePhysics> libMesh::DifferentiablePhysics::clone_physics ( )
pure virtualinherited

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  
)
inlinevirtualinherited

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 123 of file diff_physics.h.

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

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

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 105 of file diff_physics.h.

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

106  {
107  return request_jacobian;
108  }
virtual bool libMesh::FEMPhysics::eulerian_residual ( bool  request_jacobian,
DiffContext context 
)
virtual

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

This function assumes that the user's time derivative equations (except for any equations involving unknown mesh xyz coordinates themselves) are expressed in an Eulerian frame of reference, and that the user is satisfied with an unstabilized convection term. Lagrangian equations will probably require overriding eulerian_residual() with a blank function; ALE or stabilized formulations will require reimplementing eulerian_residual() entirely.

Reimplemented from libMesh::DifferentiablePhysics.

const System * libMesh::DifferentiablePhysics::get_mesh_system ( ) const
inlineinherited

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 414 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_sys.

Referenced by libMesh::FEMSystem::build_context().

415 {
416  return _mesh_sys;
417 }
System * libMesh::DifferentiablePhysics::get_mesh_system ( )
inlineinherited

Returns a reference to the system with variables corresponding to mesh nodal coordinates, or NULL if the mesh is fixed.

Definition at line 420 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_sys.

421 {
422  return _mesh_sys;
423 }
unsigned int libMesh::DifferentiablePhysics::get_mesh_x_var ( ) const
inlineinherited

Returns the variable number corresponding to the mesh x coordinate. Useful for ALE calculations.

Definition at line 426 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_x_var.

Referenced by libMesh::FEMSystem::build_context().

427 {
428  return _mesh_x_var;
429 }
unsigned int libMesh::DifferentiablePhysics::get_mesh_y_var ( ) const
inlineinherited

Returns the variable number corresponding to the mesh y coordinate. Useful for ALE calculations.

Definition at line 432 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_y_var.

Referenced by libMesh::FEMSystem::build_context().

433 {
434  return _mesh_y_var;
435 }
unsigned int libMesh::DifferentiablePhysics::get_mesh_z_var ( ) const
inlineinherited

Returns the variable number corresponding to the mesh z coordinate. Useful for ALE calculations.

Definition at line 438 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_z_var.

Referenced by libMesh::FEMSystem::build_context().

439 {
440  return _mesh_z_var;
441 }
virtual void libMesh::DifferentiablePhysics::init_context ( DiffContext )
inlinevirtualinherited

Reimplemented in libMesh::FEMSystem.

Definition at line 257 of file diff_physics.h.

257 {}
virtual void libMesh::DifferentiablePhysics::init_physics ( const System sys)
virtualinherited

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
inlineinherited

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 201 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_time_evolving.

Referenced by libMesh::FEMSystem::init_context().

201  {
202  return _time_evolving[var];
203  }
virtual bool libMesh::FEMPhysics::mass_residual ( bool  request_jacobian,
DiffContext  
)
virtual

Adds a mass vector 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.

Most problems can use the reimplementation in FEMPhysics::mass_residual; few users will need to reimplement this themselves.

Reimplemented from libMesh::DifferentiablePhysics.

void libMesh::DifferentiablePhysics::set_mesh_system ( System sys)
inlinevirtualinherited

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 370 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_sys.

371 {
372  // For now we assume that we're doing fully coupled mesh motion
373 // if (sys && sys != this)
374 // libmesh_not_implemented();
375 
376  // For the foreseeable future we'll assume that we keep these
377  // Systems in the same EquationSystems
378  // libmesh_assert_equal_to (&this->get_equation_systems(),
379  // &sys->get_equation_systems());
380 
381  // And for the immediate future this code may not even work
382  libmesh_experimental();
383 
384  _mesh_sys = sys;
385 }
void libMesh::DifferentiablePhysics::set_mesh_x_var ( unsigned int  var)
inlinevirtualinherited

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 390 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_x_var.

391 {
392  _mesh_x_var = var;
393 }
void libMesh::DifferentiablePhysics::set_mesh_y_var ( unsigned int  var)
inlinevirtualinherited

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

Definition at line 398 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_y_var.

399 {
400  _mesh_y_var = var;
401 }
void libMesh::DifferentiablePhysics::set_mesh_z_var ( unsigned int  var)
inlinevirtualinherited

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

Definition at line 406 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_mesh_z_var.

407 {
408  _mesh_z_var = var;
409 }
virtual bool libMesh::DifferentiablePhysics::side_constraint ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

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 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 172 of file diff_physics.h.

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

173  {
174  return request_jacobian;
175  }
virtual bool libMesh::DifferentiablePhysics::side_mass_residual ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

Adds a mass vector 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.

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 246 of file diff_physics.h.

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

247  {
248  return request_jacobian;
249  }
virtual bool libMesh::DifferentiablePhysics::side_time_derivative ( bool  request_jacobian,
DiffContext  
)
inlinevirtualinherited

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 152 of file diff_physics.h.

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

153  {
154  return request_jacobian;
155  }
virtual void libMesh::DifferentiablePhysics::time_evolving ( unsigned int  var)
inlinevirtualinherited

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 188 of file diff_physics.h.

References libMesh::DifferentiablePhysics::_time_evolving.

188  {
189  if (_time_evolving.size() <= var)
190  _time_evolving.resize(var+1, false);
191  _time_evolving[var] = true;
192  }

Member Data Documentation

System* libMesh::DifferentiablePhysics::_mesh_sys
protectedinherited
unsigned int libMesh::DifferentiablePhysics::_mesh_x_var
protectedinherited
unsigned int libMesh::DifferentiablePhysics::_mesh_y_var
protectedinherited
unsigned int libMesh::DifferentiablePhysics::_mesh_z_var
protectedinherited
std::vector<bool> libMesh::DifferentiablePhysics::_time_evolving
protectedinherited

Stores bools to tell us which variables are evolving in time and which are just constraints

Definition at line 360 of file diff_physics.h.

Referenced by libMesh::DifferentiablePhysics::is_time_evolving(), and libMesh::DifferentiablePhysics::time_evolving().

bool libMesh::DifferentiablePhysics::compute_internal_sides
inherited

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 134 of file diff_physics.h.


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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:58:01 UTC

Hosted By:
SourceForge.net Logo