$extrastylesheet
libMesh::FEMContext Class Reference

#include <fem_context.h>

Inheritance diagram for libMesh::FEMContext:

List of all members.

Public Types

typedef const DenseSubVector
< Number > &(DiffContext::* 
diff_subsolution_getter )(unsigned int) const
typedef std::map< const
NumericVector< Number >
*, std::pair< DenseVector
< Number >, std::vector
< DenseSubVector< Number >
* > > >::iterator 
localized_vectors_iterator

Public Member Functions

 FEMContext (const System &sys)
virtual ~FEMContext ()
bool has_side_boundary_id (boundary_id_type id) const
std::vector< boundary_id_typeside_boundary_ids () const
Number interior_value (unsigned int var, unsigned int qp) const
Number side_value (unsigned int var, unsigned int qp) const
Number point_value (unsigned int var, const Point &p) const
Gradient interior_gradient (unsigned int var, unsigned int qp) const
Gradient side_gradient (unsigned int var, unsigned int qp) const
Gradient point_gradient (unsigned int var, const Point &p) const
Tensor interior_hessian (unsigned int var, unsigned int qp) const
Tensor side_hessian (unsigned int var, unsigned int qp) const
Tensor point_hessian (unsigned int var, const Point &p) const
Number fixed_interior_value (unsigned int var, unsigned int qp) const
Number fixed_side_value (unsigned int var, unsigned int qp) const
Number fixed_point_value (unsigned int var, const Point &p) const
Gradient fixed_interior_gradient (unsigned int var, unsigned int qp) const
Gradient fixed_side_gradient (unsigned int var, unsigned int qp) const
Gradient fixed_point_gradient (unsigned int var, const Point &p) const
Tensor fixed_interior_hessian (unsigned int var, unsigned int qp) const
Tensor fixed_side_hessian (unsigned int var, unsigned int qp) const
Tensor fixed_point_hessian (unsigned int var, const Point &p) const
template<typename OutputShape >
void get_element_fe (unsigned int var, FEGenericBase< OutputShape > *&fe) const
FEBaseget_element_fe (unsigned int var) const
template<typename OutputShape >
void get_element_fe (unsigned int var, FEGenericBase< OutputShape > *&fe, unsigned char dim) const
FEBaseget_element_fe (unsigned int var, unsigned char dim) const
template<typename OutputShape >
void get_side_fe (unsigned int var, FEGenericBase< OutputShape > *&fe) const
FEBaseget_side_fe (unsigned int var) const
template<typename OutputShape >
void get_side_fe (unsigned int var, FEGenericBase< OutputShape > *&fe, unsigned char dim) const
FEBaseget_side_fe (unsigned int var, unsigned char dim) const
template<typename OutputShape >
void get_edge_fe (unsigned int var, FEGenericBase< OutputShape > *&fe) const
FEBaseget_edge_fe (unsigned int var) const
template<typename OutputType >
void interior_value (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType >
void interior_values (unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_values_vector) const
template<typename OutputType >
void side_value (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType >
void side_values (unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_values_vector) const
template<typename OutputType >
void point_value (unsigned int var, const Point &p, OutputType &u) const
template<typename OutputType >
void interior_gradient (unsigned int var, unsigned int qp, OutputType &du) const
template<typename OutputType >
void interior_gradients (unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &interior_gradients_vector) const
template<typename OutputType >
void side_gradient (unsigned int var, unsigned int qp, OutputType &du) const
template<typename OutputType >
void side_gradients (unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &side_gradients_vector) const
template<typename OutputType >
void point_gradient (unsigned int var, const Point &p, OutputType &grad_u) const
template<typename OutputType >
void interior_hessian (unsigned int var, unsigned int qp, OutputType &d2u) const
template<typename OutputType >
void interior_hessians (unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
template<typename OutputType >
void side_hessian (unsigned int var, unsigned int qp, OutputType &d2u) const
template<typename OutputType >
void side_hessians (unsigned int var, const NumericVector< Number > &_system_vector, std::vector< OutputType > &d2u_vals) const
template<typename OutputType >
void point_hessian (unsigned int var, const Point &p, OutputType &hess_u) const
template<typename OutputType >
void interior_rate (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType >
void side_rate (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType >
void point_rate (unsigned int var, const Point &p, OutputType &u) const
template<typename OutputType >
void fixed_interior_value (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType >
void fixed_side_value (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType >
void fixed_point_value (unsigned int var, const Point &p, OutputType &u) const
template<typename OutputType >
void fixed_interior_gradient (unsigned int var, unsigned int qp, OutputType &grad_u) const
template<typename OutputType >
void fixed_side_gradient (unsigned int var, unsigned int qp, OutputType &grad_u) const
template<typename OutputType >
void fixed_point_gradient (unsigned int var, const Point &p, OutputType &grad_u) const
template<typename OutputType >
void fixed_interior_hessian (unsigned int var, unsigned int qp, OutputType &hess_u) const
template<typename OutputType >
void fixed_side_hessian (unsigned int var, unsigned int qp, OutputType &hess_u) const
template<typename OutputType >
void fixed_point_hessian (unsigned int var, const Point &p, OutputType &hess_u) const
template<typename OutputType >
void interior_curl (unsigned int var, unsigned int qp, OutputType &curl_u) const
template<typename OutputType >
void point_curl (unsigned int var, const Point &p, OutputType &curl_u) const
template<typename OutputType >
void interior_div (unsigned int var, unsigned int qp, OutputType &div_u) const
virtual void elem_reinit (Real theta)
virtual void elem_side_reinit (Real theta)
virtual void elem_edge_reinit (Real theta)
virtual void nonlocal_reinit (Real theta)
virtual void pre_fe_reinit (const System &, const Elem *e)
void elem_fe_reinit ()
virtual void side_fe_reinit ()
void edge_fe_reinit ()
const QBaseget_element_qrule () const
const QBaseget_side_qrule () const
const QBaseget_element_qrule (unsigned char dim) const
const QBaseget_side_qrule (unsigned char dim) const
const QBaseget_edge_qrule () const
virtual void set_mesh_system (System *sys)
const Systemget_mesh_system () const
Systemget_mesh_system ()
unsigned int get_mesh_x_var () const
void set_mesh_x_var (unsigned int x_var)
unsigned int get_mesh_y_var () const
void set_mesh_y_var (unsigned int y_var)
unsigned int get_mesh_z_var () const
void set_mesh_z_var (unsigned int z_var)
bool has_elem () const
const Elemget_elem () const
Elemget_elem ()
unsigned char get_side () const
unsigned char get_edge () const
unsigned char get_dim () const
unsigned char get_elem_dim () const
void elem_position_set (Real theta)
void elem_position_get ()
unsigned int n_vars () const
const Systemget_system () const
const DenseVector< Number > & get_elem_solution () const
DenseVector< Number > & get_elem_solution ()
const DenseSubVector< Number > & get_elem_solution (unsigned int var) const
DenseSubVector< Number > & get_elem_solution (unsigned int var)
const DenseVector< Number > & get_elem_solution_rate () const
DenseVector< Number > & get_elem_solution_rate ()
const DenseSubVector< Number > & get_elem_solution_rate (unsigned int var) const
DenseSubVector< Number > & get_elem_solution_rate (unsigned int var)
const DenseVector< Number > & get_elem_fixed_solution () const
DenseVector< Number > & get_elem_fixed_solution ()
const DenseSubVector< Number > & get_elem_fixed_solution (unsigned int var) const
DenseSubVector< Number > & get_elem_fixed_solution (unsigned int var)
const DenseVector< Number > & get_elem_residual () const
DenseVector< Number > & get_elem_residual ()
const DenseSubVector< Number > & get_elem_residual (unsigned int var) const
DenseSubVector< Number > & get_elem_residual (unsigned int var)
const DenseMatrix< Number > & get_elem_jacobian () const
DenseMatrix< Number > & get_elem_jacobian ()
const DenseSubMatrix< Number > & get_elem_jacobian (unsigned int var1, unsigned int var2) const
DenseSubMatrix< Number > & get_elem_jacobian (unsigned int var1, unsigned int var2)
const std::vector< Number > & get_qois () const
std::vector< Number > & get_qois ()
const std::vector< DenseVector
< Number > > & 
get_qoi_derivatives () const
std::vector< DenseVector
< Number > > & 
get_qoi_derivatives ()
const DenseSubVector< Number > & get_qoi_derivatives (unsigned int qoi, unsigned int var) const
DenseSubVector< Number > & get_qoi_derivatives (unsigned int qoi, unsigned int var)
const std::vector< dof_id_type > & get_dof_indices () const
std::vector< dof_id_type > & get_dof_indices ()
const std::vector< dof_id_type > & get_dof_indices (unsigned int var) const
std::vector< dof_id_type > & get_dof_indices (unsigned int var)
Real get_system_time () const
Real get_time () const
void set_time (Real time_in)
Real get_elem_solution_derivative () const
Real get_elem_solution_rate_derivative () const
Real get_fixed_solution_derivative () const
bool is_adjoint () const
bool & is_adjoint ()
void set_deltat_pointer (Real *dt)
Real get_deltat_value ()
void add_localized_vector (NumericVector< Number > &localized_vector, const System &sys)
DenseVector< Number > & get_localized_vector (const NumericVector< Number > &localized_vector)
const DenseVector< Number > & get_localized_vector (const NumericVector< Number > &localized_vector) const
DenseSubVector< Number > & get_localized_subvector (const NumericVector< Number > &localized_vector, unsigned int var)
const DenseSubVector< Number > & get_localized_subvector (const NumericVector< Number > &localized_vector, unsigned int var) const

Public Attributes

System_mesh_sys
unsigned int _mesh_x_var
unsigned int _mesh_y_var
unsigned int _mesh_z_var
unsigned char side
unsigned char edge
Real time
const Real system_time
Real elem_solution_derivative
Real elem_solution_rate_derivative
Real fixed_solution_derivative

Protected Member Functions

template<typename OutputShape >
UniquePtr< FEGenericBase
< OutputShape > > 
build_new_fe (const FEGenericBase< OutputShape > *fe, const Point &p) const
void set_elem (const Elem *e)
template<typename OutputType , diff_subsolution_getter subsolution_getter>
void some_interior_value (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType , diff_subsolution_getter subsolution_getter>
void some_interior_gradient (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType , diff_subsolution_getter subsolution_getter>
void some_interior_hessian (unsigned int var, unsigned int qp, OutputType &u) const
template<typename OutputType , diff_subsolution_getter subsolution_getter>
void some_side_value (unsigned int var, unsigned int qp, OutputType &u) const

Protected Attributes

std::vector< std::map< FEType,
FEAbstract * > > 
_element_fe
std::vector< std::map< FEType,
FEAbstract * > > 
_side_fe
std::map< FEType, FEAbstract * > _edge_fe
std::vector< std::vector
< FEAbstract * > > 
_element_fe_var
std::vector< std::vector
< FEAbstract * > > 
_side_fe_var
std::vector< FEAbstract * > _edge_fe_var
const BoundaryInfo_boundary_info
const Elem_elem
unsigned char _dim
unsigned char _elem_dim
std::vector< QBase * > _element_qrule
std::vector< QBase * > _side_qrule
QBase_edge_qrule
std::map< const NumericVector
< Number > *, std::pair
< DenseVector< Number >
, std::vector< DenseSubVector
< Number > * > > > 
_localized_vectors
DenseVector< Number_elem_solution
std::vector< DenseSubVector
< Number > * > 
_elem_subsolutions
DenseVector< Number_elem_solution_rate
std::vector< DenseSubVector
< Number > * > 
_elem_subsolution_rates
DenseVector< Number_elem_fixed_solution
std::vector< DenseSubVector
< Number > * > 
_elem_fixed_subsolutions
DenseVector< Number_elem_residual
DenseMatrix< Number_elem_jacobian
std::vector< Number_elem_qoi
std::vector< DenseVector
< Number > > 
_elem_qoi_derivative
std::vector< std::vector
< DenseSubVector< Number > * > > 
_elem_qoi_subderivatives
std::vector< DenseSubVector
< Number > * > 
_elem_subresiduals
std::vector< std::vector
< DenseSubMatrix< Number > * > > 
_elem_subjacobians
std::vector< dof_id_type_dof_indices
std::vector< std::vector
< dof_id_type > > 
_dof_indices_var

Private Member Functions

void _do_elem_position_set (Real theta)
void _update_time_from_system (Real theta)

Detailed Description

This class provides all data required for a physics package (e.g. an FEMSystem subclass) to perform local element residual and jacobian integrations.

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 2009

Definition at line 64 of file fem_context.h.


Member Typedef Documentation

typedef const DenseSubVector<Number>&(DiffContext::* libMesh::FEMContext::diff_subsolution_getter)(unsigned int) const

Helper typedef to simplify refactoring

Definition at line 801 of file fem_context.h.

typedef std::map<const NumericVector<Number>*, std::pair<DenseVector<Number>, std::vector<DenseSubVector<Number>*> > >::iterator libMesh::DiffContext::localized_vectors_iterator [inherited]

Typedef for the localized_vectors iterator

Definition at line 480 of file diff_context.h.


Constructor & Destructor Documentation

libMesh::FEMContext::FEMContext ( const System sys) [explicit]

Constructor. Allocates some but fills no data structures.

Definition at line 36 of file fem_context.C.

References _dim, _edge_fe, _edge_fe_var, _edge_qrule, _element_fe, _element_fe_var, _element_qrule, _side_fe, _side_fe_var, _side_qrule, libMesh::FEAbstract::build(), libMesh::FEType::default_quadrature_rule(), libMesh::MeshBase::elem_dimensions(), libMesh::System::extra_quadrature_order, libMesh::FEType::family, libMesh::System::get_mesh(), libMesh::libmesh_assert(), libMesh::System::n_vars(), libMesh::FEType::order, libMesh::SCALAR, and libMesh::System::variable_type().

  : DiffContext(sys),
    _mesh_sys(NULL),
    _mesh_x_var(0),
    _mesh_y_var(0),
    _mesh_z_var(0),
    side(0), edge(0),
    _boundary_info(sys.get_mesh().get_boundary_info()),
    _elem(NULL),
    _dim(sys.get_mesh().mesh_dimension()),
    _elem_dim(0), /* This will be reset in set_elem(). */
    _element_qrule(4,NULL),
    _side_qrule(4,NULL),
    _edge_qrule(NULL)
{
  // Reserve space for the FEAbstract and QBase objects for each
  // element dimension possiblity (0,1,2,3)
  _element_fe.resize(4);
  _side_fe.resize(4);
  _element_fe_var.resize(4);
  _side_fe_var.resize(4);

  // We need to know which of our variables has the hardest
  // shape functions to numerically integrate.

  unsigned int nv = sys.n_vars();

  libmesh_assert (nv);
  FEType hardest_fe_type = sys.variable_type(0);

  bool have_scalar = false;

  for (unsigned int i=0; i != nv; ++i)
    {
      FEType fe_type = sys.variable_type(i);

      // FIXME - we don't yet handle mixed finite elements from
      // different families which require different quadrature rules
      // libmesh_assert_equal_to (fe_type.family, hardest_fe_type.family);

      if (fe_type.order > hardest_fe_type.order)
        hardest_fe_type = fe_type;

      // We need to detect SCALAR's so we can prepare FE objects for them.
      if( fe_type.family == SCALAR )
        have_scalar = true;
    }

  std::set<unsigned char> elem_dims( sys.get_mesh().elem_dimensions() );
  if(have_scalar)
    // SCALAR FEs have dimension 0 by assumption
    elem_dims.insert(0);

  for( std::set<unsigned char>::const_iterator dim_it = elem_dims.begin();
       dim_it != elem_dims.end(); ++dim_it )
    {
      const unsigned char dim = *dim_it;

      // Create an adequate quadrature rule
      _element_qrule[dim] = hardest_fe_type.default_quadrature_rule
        (dim, sys.extra_quadrature_order).release();
      _side_qrule[dim] = hardest_fe_type.default_quadrature_rule
        (dim-1, sys.extra_quadrature_order).release();
      if (dim == 3)
        _edge_qrule = hardest_fe_type.default_quadrature_rule
          (1, sys.extra_quadrature_order).release();

      // Next, create finite element objects
      _element_fe_var[dim].resize(nv);
      _side_fe_var[dim].resize(nv);
      if (dim == 3)
        _edge_fe_var.resize(nv);


      for (unsigned int i=0; i != nv; ++i)
        {
          FEType fe_type = sys.variable_type(i);

          if ( _element_fe[dim][fe_type] == NULL )
            {
              _element_fe[dim][fe_type] = FEAbstract::build(dim, fe_type).release();
              _element_fe[dim][fe_type]->attach_quadrature_rule(_element_qrule[dim]);
              _side_fe[dim][fe_type] = FEAbstract::build(dim, fe_type).release();
              _side_fe[dim][fe_type]->attach_quadrature_rule(_side_qrule[dim]);

              if (this->_dim == 3)
                {
                  _edge_fe[fe_type] = FEAbstract::build(dim, fe_type).release();
                  _edge_fe[fe_type]->attach_quadrature_rule(_edge_qrule);
                }
            }

          _element_fe_var[dim][i] = _element_fe[dim][fe_type];
          _side_fe_var[dim][i] = _side_fe[dim][fe_type];
          if ((dim) == 3)
            _edge_fe_var[i] = _edge_fe[fe_type];

        }
    }
}

Destructor.

Definition at line 137 of file fem_context.C.

References _edge_fe, _edge_qrule, _element_fe, _element_qrule, _side_fe, and _side_qrule.

{
  // We don't want to store UniquePtrs in STL containers, but we don't
  // want to leak memory either
  for (std::vector<std::map<FEType, FEAbstract *> >::iterator d = _element_fe.begin();
       d != _element_fe.end(); ++d)
    for (std::map<FEType, FEAbstract *>::iterator i = d->begin();
         i != d->end(); ++i)
      delete i->second;

  for (std::vector<std::map<FEType, FEAbstract *> >::iterator d = _side_fe.begin();
       d != _side_fe.end(); ++d)
    for (std::map<FEType, FEAbstract *>::iterator i = d->begin();
         i != d->end(); ++i)
      delete i->second;

  for (std::map<FEType, FEAbstract *>::iterator i = _edge_fe.begin();
       i != _edge_fe.end(); ++i)
    delete i->second;
  _edge_fe.clear();

  for (std::vector<QBase*>::iterator i = _element_qrule.begin();
       i != _element_qrule.end(); ++i)
    delete *i;
  _element_qrule.clear();

  for (std::vector<QBase*>::iterator i = _side_qrule.begin();
       i != _side_qrule.end(); ++i)
    delete *i;
  _side_qrule.clear();

  delete _edge_qrule;
  _edge_qrule = NULL;
}

Member Function Documentation

Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to the value it would take after a fraction theta of a timestep.

This does the work of elem_position_set, but isn't safe to call without _mesh_sys/etc. defined first.

Definition at line 1557 of file fem_context.C.

References _mesh_sys, get_elem(), get_elem_dim(), libMesh::DiffContext::get_elem_solution(), get_element_fe(), get_mesh_x_var(), get_mesh_y_var(), get_mesh_z_var(), libMesh::invalid_uint, libMesh::LAGRANGE, libMesh::libmesh_assert(), libMesh::libmesh_real(), n_nodes, libMesh::Elem::n_nodes(), and libMesh::n_threads().

Referenced by elem_position_set().

{
  // This is too expensive to call unless we've been asked to move the mesh
  libmesh_assert (_mesh_sys);

  // This will probably break with threading when two contexts are
  // operating on elements which share a node
  libmesh_assert_equal_to (libMesh::n_threads(), 1);

  const unsigned char dim = this->get_elem_dim();

  // If the coordinate data is in our own system, it's already
  // been set up for us, and we can ignore our input parameter theta
  //  if (_mesh_sys == this->number())
  //    {
  unsigned int n_nodes = this->get_elem().n_nodes();
  // For simplicity we demand that mesh coordinates be stored
  // in a format that allows a direct copy
  libmesh_assert(this->get_mesh_x_var() == libMesh::invalid_uint ||
                 (this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().family
                  == LAGRANGE &&
                  this->get_elem_solution(this->get_mesh_x_var()).size() == n_nodes));
  libmesh_assert(this->get_mesh_y_var() == libMesh::invalid_uint ||
                 (this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().family
                  == LAGRANGE &&
                  this->get_elem_solution(this->get_mesh_y_var()).size() == n_nodes));
  libmesh_assert(this->get_mesh_z_var() == libMesh::invalid_uint ||
                 (this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().family
                  == LAGRANGE &&
                  this->get_elem_solution(this->get_mesh_z_var()).size() == n_nodes));

  // Set the new point coordinates
  if (this->get_mesh_x_var() != libMesh::invalid_uint)
    for (unsigned int i=0; i != n_nodes; ++i)
      const_cast<Elem&>(this->get_elem()).point(i)(0) =
        libmesh_real(this->get_elem_solution(this->get_mesh_x_var())(i));

  if (this->get_mesh_y_var() != libMesh::invalid_uint)
    for (unsigned int i=0; i != n_nodes; ++i)
      const_cast<Elem&>(this->get_elem()).point(i)(1) =
        libmesh_real(this->get_elem_solution(this->get_mesh_y_var())(i));

  if (this->get_mesh_z_var() != libMesh::invalid_uint)
    for (unsigned int i=0; i != n_nodes; ++i)
      const_cast<Elem&>(this->get_elem()).point(i)(2) =
        libmesh_real(this->get_elem_solution(this->get_mesh_z_var())(i));
  //    }
  // FIXME - If the coordinate data is not in our own system, someone
  // had better get around to implementing that... - RHS
  //  else
  //    {
  //      libmesh_not_implemented();
  //    }
}

Update the time in the context object for the given value of theta, based on the values of "time" and "deltat" stored in the system which created this context.

Definition at line 1751 of file fem_context.C.

References libMesh::DiffContext::get_deltat_value(), libMesh::DiffContext::get_system_time(), libMesh::Real, and libMesh::DiffContext::set_time().

Referenced by elem_edge_reinit(), elem_reinit(), elem_side_reinit(), and nonlocal_reinit().

{
  // Update the "time" variable based on the value of theta.  For this
  // to work, we need to know the value of deltat, a pointer to which is now
  // stored by our parent DiffContext class.  Note: get_deltat_value() will
  // assert in debug mode if the requested pointer is NULL.
  const Real deltat = this->get_deltat_value();

  this->set_time(theta*(this->get_system_time() + deltat) + (1.-theta)*this->get_system_time());
}
void libMesh::DiffContext::add_localized_vector ( NumericVector< Number > &  localized_vector,
const System sys 
) [inherited]

Adds a vector to the map of localized vectors. We can later evaluate interior_values, interior_gradients and side_values for these fields these vectors represent.

Definition at line 129 of file diff_context.C.

References libMesh::DiffContext::_localized_vectors, and libMesh::System::n_vars().

{
  // Make an empty pair keyed with a reference to this _localized_vector
  _localized_vectors[&localized_vector] = std::make_pair(DenseVector<Number>(), std::vector<DenseSubVector<Number>*>());

  unsigned int nv = sys.n_vars();

  _localized_vectors[&localized_vector].second.reserve(nv);

  // Fill the DenseSubVector with nv copies of DenseVector
  for(unsigned int i=0; i != nv; ++i)
    _localized_vectors[&localized_vector].second.push_back(new DenseSubVector<Number>(_localized_vectors[&localized_vector].first));
}
template<typename OutputShape >
template UniquePtr< FEGenericBase< RealGradient > > libMesh::FEMContext::build_new_fe ( const FEGenericBase< OutputShape > *  fe,
const Point p 
) const [protected]

Helper function to reduce some code duplication in the *_point_* methods.

Definition at line 1765 of file fem_context.C.

References libMesh::Elem::dim(), libMesh::FEType::family, get_elem(), libMesh::FEAbstract::get_fe_type(), has_elem(), libMesh::FEInterface::inverse_map(), libMesh::libmesh_assert(), libMesh::FEAbstract::reinit(), and libMesh::SCALAR.

Referenced by fixed_point_gradient(), fixed_point_hessian(), fixed_point_value(), point_curl(), point_gradient(), point_hessian(), and point_value().

{
  FEType fe_type = fe->get_fe_type();

  // If we don't have an Elem to evaluate on, then the only functions
  // we can sensibly evaluate are the scalar dofs which are the same
  // everywhere.
  libmesh_assert(this->has_elem() || fe_type.family == SCALAR);

  unsigned int elem_dim = this->has_elem() ? this->get_elem().dim() : 0;

  FEGenericBase<OutputShape>* fe_new =
    FEGenericBase<OutputShape>::build(elem_dim, fe_type).release();

  // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h
  // Build a vector of point co-ordinates to send to reinit
  Point master_point = this->has_elem() ?
    FEInterface::inverse_map(elem_dim, fe_type, &this->get_elem(), p) :
    Point(0);

  std::vector<Point> coor(1, master_point);

  // Reinitialize the element and compute the shape function values at coor
  if(this->has_elem())
    fe_new->reinit (&this->get_elem(), &coor);
  else
    // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
    fe_new->reinit (NULL, &coor);

  return UniquePtr<FEGenericBase<OutputShape> >(fe_new);
}

Reinitializes edge FE objects on the current geometric element

Definition at line 1478 of file fem_context.C.

References _edge_fe, get_edge(), get_elem(), and get_elem_dim().

Referenced by elem_edge_reinit().

{
  libmesh_assert_equal_to (this->get_elem_dim(), 3);

  // Initialize all the interior FE objects on elem/edge.
  // Logging of FE::reinit is done in the FE functions
  std::map<FEType, FEAbstract *>::iterator local_fe_end = _edge_fe.end();
  for (std::map<FEType, FEAbstract *>::iterator i = _edge_fe.begin();
       i != local_fe_end; ++i)
    {
      i->second->edge_reinit(&(this->get_elem()), this->get_edge());
    }
}
void libMesh::FEMContext::elem_edge_reinit ( Real  theta) [virtual]

Reinitialize Elem and edge FE objects if necessary for integration at a new point in time: specifically, handle moving elements in moving mesh schemes.

Reimplemented from libMesh::DiffContext.

Definition at line 1410 of file fem_context.C.

References _mesh_sys, _update_time_from_system(), edge_fe_reinit(), and elem_position_set().

{
  // Update the "time" variable of this context object
  this->_update_time_from_system(theta);

  // Handle a moving element if necessary
  if (_mesh_sys)
    {
      // FIXME - not threadsafe yet!
      elem_position_set(theta);
      edge_fe_reinit();
    }
}

Reinitializes interior FE objects on the current geometric element

Definition at line 1435 of file fem_context.C.

References _element_fe, get_elem(), get_elem_dim(), has_elem(), and libMesh::libmesh_assert().

Referenced by elem_reinit(), libMesh::FEMSystem::mesh_position_set(), and nonlocal_reinit().

{
  // Initialize all the interior FE objects on elem.
  // Logging of FE::reinit is done in the FE functions
  // We only reinit the FE objects for the current element
  // dimension
  const unsigned char dim = this->get_elem_dim();

  libmesh_assert( !_element_fe[dim].empty() );

  std::map<FEType, FEAbstract *>::iterator local_fe_end = _element_fe[dim].end();
  for (std::map<FEType, FEAbstract *>::iterator i = _element_fe[dim].begin();
       i != local_fe_end; ++i)
    {
      if(this->has_elem())
        i->second->reinit(&(this->get_elem()));
      else
        // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
        i->second->reinit(NULL);
    }
}

Uses the geometry of elem to set the coordinate data specified by mesh_*_position configuration.

Definition at line 1494 of file fem_context.C.

References _mesh_sys, get_elem(), get_elem_dim(), libMesh::DiffContext::get_elem_solution(), get_element_fe(), get_mesh_x_var(), get_mesh_y_var(), get_mesh_z_var(), libMesh::invalid_uint, libMesh::LAGRANGE, libMesh::libmesh_assert(), n_nodes, libMesh::Elem::n_nodes(), libMesh::n_threads(), and libMesh::Elem::point().

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

{
  // This is too expensive to call unless we've been asked to move the mesh
  libmesh_assert (_mesh_sys);

  // This will probably break with threading when two contexts are
  // operating on elements which share a node
  libmesh_assert_equal_to (libMesh::n_threads(), 1);

  // If the coordinate data is in our own system, it's already
  // been set up for us
  //  if (_mesh_sys == this->number())
  //    {
  unsigned int n_nodes = this->get_elem().n_nodes();

  const unsigned char dim = this->get_elem_dim();

  // For simplicity we demand that mesh coordinates be stored
  // in a format that allows a direct copy
  libmesh_assert(this->get_mesh_x_var() == libMesh::invalid_uint ||
                 (this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().family
                  == LAGRANGE &&
                  this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().order
                  == this->get_elem().default_order()));
  libmesh_assert(this->get_mesh_y_var() == libMesh::invalid_uint ||
                 (this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().family
                  == LAGRANGE &&
                  this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().order
                  == this->get_elem().default_order()));
  libmesh_assert(this->get_mesh_z_var() == libMesh::invalid_uint ||
                 (this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().family
                  == LAGRANGE &&
                  this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().order
                  == this->get_elem().default_order()));

  // Get degree of freedom coefficients from point coordinates
  if (this->get_mesh_x_var() != libMesh::invalid_uint)
    for (unsigned int i=0; i != n_nodes; ++i)
      (this->get_elem_solution(this->get_mesh_x_var()))(i) = this->get_elem().point(i)(0);

  if (this->get_mesh_y_var() != libMesh::invalid_uint)
    for (unsigned int i=0; i != n_nodes; ++i)
      (this->get_elem_solution(this->get_mesh_y_var()))(i) = this->get_elem().point(i)(1);

  if (this->get_mesh_z_var() != libMesh::invalid_uint)
    for (unsigned int i=0; i != n_nodes; ++i)
      (this->get_elem_solution(this->get_mesh_z_var()))(i) = this->get_elem().point(i)(2);
  //    }
  // FIXME - If the coordinate data is not in our own system, someone
  // had better get around to implementing that... - RHS
  //  else
  //    {
  //      libmesh_not_implemented();
  //    }
}
void libMesh::FEMContext::elem_position_set ( Real  theta) [inline]

Uses the coordinate data specified by mesh_*_position configuration to set the geometry of elem to the value it would take after a fraction theta of a timestep.

Definition at line 930 of file fem_context.h.

References _do_elem_position_set(), and _mesh_sys.

Referenced by elem_edge_reinit(), elem_reinit(), elem_side_reinit(), and libMesh::FEMSystem::mesh_position_set().

{
  if (_mesh_sys)
    this->_do_elem_position_set(theta);
}
void libMesh::FEMContext::elem_reinit ( Real  theta) [virtual]

Reinitialize all my FEM context data on a given element for the given system Reinitialize Elem and FE objects if necessary for integration at a new point in time: specifically, handle moving elements in moving mesh schemes.

Reimplemented from libMesh::DiffContext.

Definition at line 1373 of file fem_context.C.

References _mesh_sys, _update_time_from_system(), elem_fe_reinit(), elem_position_set(), and libMesh::n_threads().

{
  // Update the "time" variable of this context object
  this->_update_time_from_system(theta);

  // Handle a moving element if necessary.
  if (_mesh_sys)
    // We assume that the ``default'' state
    // of the mesh is its final, theta=1.0
    // position, so we don't bother with
    // mesh motion in that case.
    if (theta != 1.0)
      {
        // FIXME - ALE is not threadsafe yet!
        libmesh_assert_equal_to (libMesh::n_threads(), 1);

        elem_position_set(theta);
      }
  elem_fe_reinit();
}
void libMesh::FEMContext::elem_side_reinit ( Real  theta) [virtual]

Reinitialize Elem and side FE objects if necessary for integration at a new point in time: specifically, handle moving elements in moving mesh schemes.

Reimplemented from libMesh::DiffContext.

Definition at line 1395 of file fem_context.C.

References _mesh_sys, _update_time_from_system(), elem_position_set(), and side_fe_reinit().

{
  // Update the "time" variable of this context object
  this->_update_time_from_system(theta);

  // Handle a moving element if necessary
  if (_mesh_sys)
    {
      // FIXME - not threadsafe yet!
      elem_position_set(theta);
      side_fe_reinit();
    }
}
Gradient libMesh::FEMContext::fixed_interior_gradient ( unsigned int  var,
unsigned int  qp 
) const

Returns the gradient of the fixed_solution variable var at the quadrature point qp on the current element interior. This API currently present for backward compatibility.

Definition at line 1034 of file fem_context.C.

{
  Gradient du;

  this->fixed_interior_gradient( var, qp, du );

  return du;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_interior_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  grad_u 
) const

Returns the gradient of the fixed_solution variable var at the quadrature point qp on the current element interior. This is the preferred API.

Tensor libMesh::FEMContext::fixed_interior_hessian ( unsigned int  var,
unsigned int  qp 
) const

Returns the hessian of the fixed_solution variable var at the quadrature point qp on the current element interior. This API currently present for backward compatibility.

Definition at line 1055 of file fem_context.C.

{
  Tensor d2u;

  this->fixed_interior_hessian( var, qp, d2u );

  return d2u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_interior_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  hess_u 
) const

Returns the hessian of the fixed_solution variable var at the quadrature point qp on the current element interior. This is the preferred API.

Definition at line 1066 of file fem_context.C.

{
  this->some_interior_hessian <OutputType,&DiffContext::get_elem_fixed_solution>
    (var, qp, d2u);
}
Number libMesh::FEMContext::fixed_interior_value ( unsigned int  var,
unsigned int  qp 
) const

Returns the value of the fixed_solution variable var at the quadrature point qp on the current element interior. This API currently present for backward compatibility.

Definition at line 1013 of file fem_context.C.

{
  Number u = 0.;

  this->fixed_interior_value( var, qp, u );

  return u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_interior_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const

Returns the value of the fixed_solution variable var at the quadrature point qp on the current element interior. This is the preferred API.

Definition at line 1025 of file fem_context.C.

{
  this->some_interior_value <OutputType,&DiffContext::get_elem_fixed_solution>
    (var, qp, u);
}
Gradient libMesh::FEMContext::fixed_point_gradient ( unsigned int  var,
const Point p 
) const

Returns the gradient of the fixed_solution variable var at the physical point p on the current element. This API currently present for backward compatibility.

Definition at line 1260 of file fem_context.C.

{
  Gradient grad_u;

  this->fixed_point_gradient( var, p, grad_u );

  return grad_u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_point_gradient ( unsigned int  var,
const Point p,
OutputType &  grad_u 
) const

Returns the gradient of the fixed_solution variable var at the physical point p on the current element. This is the preferred API.

Definition at line 1272 of file fem_context.C.

References libMesh::DiffContext::_elem_fixed_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_fixed_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >&  dphi = fe_new->get_dphi();

  grad_u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    grad_u.add_scaled(dphi[l][0], coef(l));

  return;
}
Tensor libMesh::FEMContext::fixed_point_hessian ( unsigned int  var,
const Point p 
) const

Returns the hessian of the fixed_solution variable var at the physical point p on the current element. This API currently present for backward compatibility.

Definition at line 1310 of file fem_context.C.

{
  Tensor hess_u;

  this->fixed_point_hessian( var, p, hess_u );

  return hess_u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_point_hessian ( unsigned int  var,
const Point p,
OutputType &  hess_u 
) const

Returns the hessian of the fixed_solution variable var at the physical point p on the current element. This is the preferred API.

Definition at line 1322 of file fem_context.C.

References libMesh::DiffContext::_elem_fixed_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_fixed_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >&  d2phi = fe_new->get_d2phi();

  hess_u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    hess_u.add_scaled(d2phi[l][0], coef(l));

  return;
}
Number libMesh::FEMContext::fixed_point_value ( unsigned int  var,
const Point p 
) const

Returns the value of the fixed_solution variable var at the physical point p on the current element. This API currently present for backward compatibility.

Definition at line 1215 of file fem_context.C.

{
  Number u = 0.;

  this->fixed_point_value( var, p, u );

  return u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_point_value ( unsigned int  var,
const Point p,
OutputType &  u 
) const

Returns the value of the fixed_solution variable var at the physical point p on the current element. This is the preferred API.

Definition at line 1225 of file fem_context.C.

References libMesh::DiffContext::_elem_fixed_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_fixed_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<OutputShape> >&  phi = fe_new->get_phi();

  u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    u += phi[l][0] * coef(l);

  return;
}
Gradient libMesh::FEMContext::fixed_side_gradient ( unsigned int  var,
unsigned int  qp 
) const

Returns the gradient of the fixed_solution variable var at the quadrature point qp on the current element side. This API currently present for backward compatibility.

Definition at line 1120 of file fem_context.C.

{
  Gradient du;

  this->fixed_side_gradient( var, qp, du );

  return du;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_side_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  grad_u 
) const

Returns the gradient of the fixed_solution variable var at the quadrature point qp on the current element side. This is the preferred API.

Definition at line 1131 of file fem_context.C.

References libMesh::DiffContext::_elem_fixed_subsolutions, libMesh::DiffContext::get_dof_indices(), libMesh::FEGenericBase< OutputType >::get_dphi(), libMesh::DiffContext::get_elem_fixed_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_fixed_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();

  // Accumulate solution derivatives
  du = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    du.add_scaled(dphi[l][qp], coef(l));

  return;
}
Tensor libMesh::FEMContext::fixed_side_hessian ( unsigned int  var,
unsigned int  qp 
) const

Returns the hessian of the fixed_solution variable var at the quadrature point qp on the current element side. This API currently present for backward compatibility.

Definition at line 1167 of file fem_context.C.

{
  Tensor d2u;

  this->fixed_side_hessian( var, qp, d2u );

  return d2u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_side_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  hess_u 
) const

Returns the hessian of the fixed_solution variable var at the quadrature point qp on the current element side. This is the preferred API.

Definition at line 1177 of file fem_context.C.

References libMesh::DiffContext::_elem_fixed_subsolutions, libMesh::FEGenericBase< OutputType >::get_d2phi(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_fixed_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();

  // Accumulate solution second derivatives
  d2u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    d2u.add_scaled(d2phi[l][qp], coef(l));

  return;
}
Number libMesh::FEMContext::fixed_side_value ( unsigned int  var,
unsigned int  qp 
) const

Returns the value of the fixed_solution variable var at the quadrature point qp on the current element side. This API currently present for backward compatibility.

Definition at line 1076 of file fem_context.C.

{
  Number u = 0.;

  this->fixed_side_value( var, qp, u );

  return u;
}
template<typename OutputType >
void libMesh::FEMContext::fixed_side_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const

Returns the value of the fixed_solution variable var at the quadrature point qp on the current element side. This is the preferred API.

Definition at line 1087 of file fem_context.C.

References libMesh::DiffContext::_elem_fixed_subsolutions, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), libMesh::FEGenericBase< OutputType >::get_phi(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_fixed_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();

  // Accumulate solution value
  u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    u += phi[l][qp] * coef(l);

  return;
}

Returns the value currently pointed to by this class's _deltat member

Definition at line 121 of file diff_context.C.

References libMesh::DiffContext::_deltat, and libMesh::libmesh_assert().

Referenced by _update_time_from_system().

unsigned char libMesh::FEMContext::get_dim ( ) const [inline]

Accessor for cached mesh dimension. This is the largest dimension of the elements in the mesh. For the dimension of this->_elem, use get_elem_dim();

Definition at line 737 of file fem_context.h.

References _dim.

Referenced by get_element_fe(), get_side_fe(), and libMesh::DGFEMContext::neighbor_side_fe_reinit().

  { return this->_dim; }
std::vector<dof_id_type>& libMesh::DiffContext::get_dof_indices ( ) [inline, inherited]

Non-const accessor for element dof indices

Definition at line 341 of file diff_context.h.

References libMesh::DiffContext::_dof_indices.

  { return _dof_indices; }
const std::vector<dof_id_type>& libMesh::DiffContext::get_dof_indices ( unsigned int  var) const [inline, inherited]

Accessor for element dof indices of a particular variable corresponding to the index argument.

Definition at line 348 of file diff_context.h.

References libMesh::DiffContext::_dof_indices_var.

  {
    libmesh_assert_greater(_dof_indices_var.size(), var);
    return _dof_indices_var[var];
  }
std::vector<dof_id_type>& libMesh::DiffContext::get_dof_indices ( unsigned int  var) [inline, inherited]

Accessor for element dof indices of a particular variable corresponding to the index argument.

Definition at line 358 of file diff_context.h.

References libMesh::DiffContext::_dof_indices_var.

  {
    libmesh_assert_greater(_dof_indices_var.size(), var);
    return _dof_indices_var[var];
  }
unsigned char libMesh::FEMContext::get_edge ( ) const [inline]

Accessor for current edge of Elem object

Definition at line 729 of file fem_context.h.

References edge.

Referenced by edge_fe_reinit().

  { return edge; }
template<typename OutputShape >
void libMesh::FEMContext::get_edge_fe ( unsigned int  var,
FEGenericBase< OutputShape > *&  fe 
) const [inline]

Accessor for edge (3D only!) finite element object for variable var.

Definition at line 974 of file fem_context.h.

References _edge_fe_var.

{
  libmesh_assert_less ( var, _edge_fe_var.size() );
  fe = cast_ptr<FEGenericBase<OutputShape>*>( _edge_fe_var[var] );
}
FEBase * libMesh::FEMContext::get_edge_fe ( unsigned int  var) const [inline]

Accessor for edge (3D only!) finite element object for scalar-valued variable var.

Definition at line 981 of file fem_context.h.

References _edge_fe_var.

{
  libmesh_assert_less ( var, _edge_fe_var.size() );
  return cast_ptr<FEBase*>( _edge_fe_var[var] );
}
const QBase& libMesh::FEMContext::get_edge_qrule ( ) const [inline]

Accessor for element edge quadrature rule.

Definition at line 632 of file fem_context.h.

References _edge_qrule.

  { return *(this->_edge_qrule); }

Accessor for current Elem object

Definition at line 716 of file fem_context.h.

References _elem, and libMesh::libmesh_assert().

  { libmesh_assert(this->_elem);
    return *(const_cast<Elem*>(this->_elem)); }
unsigned char libMesh::FEMContext::get_elem_dim ( ) const [inline]

Returns the dimension of this->_elem. For mixed dimension meshes, this may be different from get_dim().

Definition at line 744 of file fem_context.h.

References _elem_dim.

Referenced by _do_elem_position_set(), edge_fe_reinit(), elem_fe_reinit(), elem_position_get(), get_element_qrule(), get_side_qrule(), and side_fe_reinit().

  { return _elem_dim; }

Non-const accessor for element fixed solution.

Definition at line 189 of file diff_context.h.

References libMesh::DiffContext::_elem_fixed_solution.

  { return _elem_fixed_solution; }
const DenseSubVector<Number>& libMesh::DiffContext::get_elem_fixed_solution ( unsigned int  var) const [inline, inherited]

Accessor for element fixed solution of a particular variable corresponding to the variable index argument.

Definition at line 196 of file diff_context.h.

References libMesh::DiffContext::_elem_fixed_subsolutions, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_fixed_subsolutions.size(), var);
    libmesh_assert(_elem_fixed_subsolutions[var]);
    return *(_elem_fixed_subsolutions[var]);
  }
DenseSubVector<Number>& libMesh::DiffContext::get_elem_fixed_solution ( unsigned int  var) [inline, inherited]

Accessor for element fixed solution of a particular variable corresponding to the variable index argument.

Definition at line 207 of file diff_context.h.

References libMesh::DiffContext::_elem_fixed_subsolutions, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_fixed_subsolutions.size(), var);
    libmesh_assert(_elem_fixed_subsolutions[var]);
    return *(_elem_fixed_subsolutions[var]);
  }

Non-const accessor for element Jacobian.

Definition at line 257 of file diff_context.h.

References libMesh::DiffContext::_elem_jacobian.

  { return _elem_jacobian; }
const DenseSubMatrix<Number>& libMesh::DiffContext::get_elem_jacobian ( unsigned int  var1,
unsigned int  var2 
) const [inline, inherited]

Const accessor for element Jacobian of particular variables corresponding to the variable index arguments.

Definition at line 264 of file diff_context.h.

References libMesh::DiffContext::_elem_subjacobians, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subjacobians.size(), var1);
    libmesh_assert_greater(_elem_subjacobians[var1].size(), var2);
    libmesh_assert(_elem_subjacobians[var1][var2]);
    return *(_elem_subjacobians[var1][var2]);
  }
DenseSubMatrix<Number>& libMesh::DiffContext::get_elem_jacobian ( unsigned int  var1,
unsigned int  var2 
) [inline, inherited]

Non-const accessor for element Jacobian of particular variables corresponding to the variable index arguments.

Definition at line 276 of file diff_context.h.

References libMesh::DiffContext::_elem_subjacobians, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subjacobians.size(), var1);
    libmesh_assert_greater(_elem_subjacobians[var1].size(), var2);
    libmesh_assert(_elem_subjacobians[var1][var2]);
    return *(_elem_subjacobians[var1][var2]);
  }

Non-const accessor for element residual.

Definition at line 223 of file diff_context.h.

References libMesh::DiffContext::_elem_residual.

  { return _elem_residual; }
const DenseSubVector<Number>& libMesh::DiffContext::get_elem_residual ( unsigned int  var) const [inline, inherited]

Const accessor for element residual of a particular variable corresponding to the variable index argument.

Definition at line 230 of file diff_context.h.

References libMesh::DiffContext::_elem_subresiduals, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subresiduals.size(), var);
    libmesh_assert(_elem_subresiduals[var]);
    return *(_elem_subresiduals[var]);
  }
DenseSubVector<Number>& libMesh::DiffContext::get_elem_residual ( unsigned int  var) [inline, inherited]

Non-const accessor for element residual of a particular variable corresponding to the variable index argument.

Definition at line 241 of file diff_context.h.

References libMesh::DiffContext::_elem_subresiduals, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subresiduals.size(), var);
    libmesh_assert(_elem_subresiduals[var]);
    return *(_elem_subresiduals[var]);
  }

Non-const accessor for element solution.

Definition at line 119 of file diff_context.h.

References libMesh::DiffContext::_elem_solution.

  { return _elem_solution; }
const DenseSubVector<Number>& libMesh::DiffContext::get_elem_solution ( unsigned int  var) const [inline, inherited]

Accessor for element solution of a particular variable corresponding to the variable index argument.

Definition at line 126 of file diff_context.h.

References libMesh::DiffContext::_elem_subsolutions, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subsolutions.size(), var);
    libmesh_assert(_elem_subsolutions[var]);
    return *(_elem_subsolutions[var]);
  }
DenseSubVector<Number>& libMesh::DiffContext::get_elem_solution ( unsigned int  var) [inline, inherited]

Accessor for element solution of a particular variable corresponding to the variable index argument.

Definition at line 137 of file diff_context.h.

References libMesh::DiffContext::_elem_subsolutions, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subsolutions.size(), var);
    libmesh_assert(_elem_subsolutions[var]);
    return *(_elem_subsolutions[var]);
  }

The derivative of the current elem_solution w.r.t. the unknown solution. Corresponding Jacobian contributions should be multiplied by this amount, or may be skipped if get_elem_solution_derivative() is 0.

Definition at line 388 of file diff_context.h.

References libMesh::DiffContext::elem_solution_derivative.

const DenseVector<Number>& libMesh::DiffContext::get_elem_solution_rate ( ) const [inline, inherited]

Accessor for element solution rate of change w.r.t. time.

Definition at line 147 of file diff_context.h.

References libMesh::DiffContext::_elem_solution_rate.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), and pre_fe_reinit().

  { return _elem_solution_rate; }

Non-const accessor for element solution rate of change w.r.t. time.

Definition at line 154 of file diff_context.h.

References libMesh::DiffContext::_elem_solution_rate.

  { return _elem_solution_rate; }
const DenseSubVector<Number>& libMesh::DiffContext::get_elem_solution_rate ( unsigned int  var) const [inline, inherited]

Accessor for element solution rate for a particular variable corresponding to the variable index argument.

Definition at line 161 of file diff_context.h.

References libMesh::DiffContext::_elem_subsolution_rates, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subsolution_rates.size(), var);
    libmesh_assert(_elem_subsolution_rates[var]);
    return *(_elem_subsolution_rates[var]);
  }
DenseSubVector<Number>& libMesh::DiffContext::get_elem_solution_rate ( unsigned int  var) [inline, inherited]

Accessor for element solution rate for a particular variable corresponding to the variable index argument.

Definition at line 172 of file diff_context.h.

References libMesh::DiffContext::_elem_subsolution_rates, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_subsolution_rates.size(), var);
    libmesh_assert(_elem_subsolution_rates[var]);
    return *(_elem_subsolution_rates[var]);
  }

The derivative of the current elem_solution_rate w.r.t. the unknown solution. Corresponding Jacobian contributions should be multiplied by this amount, or may be skipped if get_elem_solution_rate_derivative() is 0.

Definition at line 397 of file diff_context.h.

References libMesh::DiffContext::elem_solution_rate_derivative.

template<typename OutputShape >
void libMesh::FEMContext::get_element_fe ( unsigned int  var,
FEGenericBase< OutputShape > *&  fe 
) const [inline]

Accessor for interior finite element object for variable var for the largest dimension in the mesh. If you have lower dimensional elements in the mesh and need to query for those FE objects, use the alternative get_element_fe method.

Definition at line 228 of file fem_context.h.

References get_dim().

Referenced by _do_elem_position_set(), elem_position_get(), libMesh::FEMSystem::init_context(), and libMesh::ParsedFEMFunction< Output >::init_context().

  { this->get_element_fe<OutputShape>(var,fe,this->get_dim()); }
FEBase* libMesh::FEMContext::get_element_fe ( unsigned int  var) const [inline]

Accessor for interior finite element object for scalar-valued variable var for the largest dimension in the mesh. If you have lower dimensional elements in the mesh and need to query for those FE objects, use the alternative get_element_fe method.

Definition at line 237 of file fem_context.h.

References get_dim(), and get_element_fe().

Referenced by get_element_fe().

  { return this->get_element_fe(var,this->get_dim()); }
template<typename OutputShape >
void libMesh::FEMContext::get_element_fe ( unsigned int  var,
FEGenericBase< OutputShape > *&  fe,
unsigned char  dim 
) const [inline]

Accessor for interior finite element object for variable var for dimension dim.

Definition at line 938 of file fem_context.h.

References _element_fe_var, and libMesh::libmesh_assert().

{
  libmesh_assert( !_element_fe_var[dim].empty() );
  libmesh_assert_less ( var, (_element_fe_var[dim].size() ) );
  fe = cast_ptr<FEGenericBase<OutputShape>*>( (_element_fe_var[dim][var] ) );
}
FEBase * libMesh::FEMContext::get_element_fe ( unsigned int  var,
unsigned char  dim 
) const [inline]

Accessor for interior finite element object for scalar-valued variable var for dimension dim.

Definition at line 947 of file fem_context.h.

References _element_fe_var, and libMesh::libmesh_assert().

{
  libmesh_assert( !_element_fe_var[dim].empty() );
  libmesh_assert_less ( var, (_element_fe_var[dim].size() ) );
  return cast_ptr<FEBase*>( (_element_fe_var[dim][var] ) );
}
const QBase& libMesh::FEMContext::get_element_qrule ( ) const [inline]

Accessor for element interior quadrature rule for the dimension of the current _elem.

Definition at line 605 of file fem_context.h.

References get_elem_dim(), and get_element_qrule().

Referenced by get_element_qrule().

  { return this->get_element_qrule(this->get_elem_dim()); }
const QBase& libMesh::FEMContext::get_element_qrule ( unsigned char  dim) const [inline]

Accessor for element interior quadrature rule.

Definition at line 618 of file fem_context.h.

References _element_qrule, and libMesh::libmesh_assert().

  { libmesh_assert(_element_qrule[dim]);
    return *(this->_element_qrule[dim]); }

The derivative of the current fixed_elem_solution w.r.t. the unknown solution. Corresponding Jacobian contributions should be multiplied by this amount, or may be skipped if get_fixed_elem_solution_derivative() is 0.

Definition at line 406 of file diff_context.h.

References libMesh::DiffContext::fixed_solution_derivative.

DenseSubVector< Number > & libMesh::DiffContext::get_localized_subvector ( const NumericVector< Number > &  localized_vector,
unsigned int  var 
) [inherited]

Return a reference to DenseSubVector localization of localized_vector at variable var contained in the _localized_vectors map

Definition at line 159 of file diff_context.C.

References libMesh::DiffContext::_localized_vectors.

Referenced by interior_values().

{
  return *_localized_vectors[&localized_vector].second[var];
}
const DenseSubVector< Number > & libMesh::DiffContext::get_localized_subvector ( const NumericVector< Number > &  localized_vector,
unsigned int  var 
) const [inherited]

const accessible version of get_localized_subvector function

Definition at line 165 of file diff_context.C.

References libMesh::DiffContext::_localized_vectors, and libMesh::libmesh_assert().

{
  std::map<const NumericVector<Number>*, std::pair<DenseVector<Number>, std::vector<DenseSubVector<Number>*> > >::const_iterator
    localized_vectors_it = _localized_vectors.find(&localized_vector);
  libmesh_assert(localized_vectors_it != _localized_vectors.end());
  return *localized_vectors_it->second.second[var];
}
DenseVector< Number > & libMesh::DiffContext::get_localized_vector ( const NumericVector< Number > &  localized_vector) [inherited]

Return a reference to DenseVector localization of localized_vector contained in the _localized_vectors map

Definition at line 144 of file diff_context.C.

References libMesh::DiffContext::_localized_vectors.

{
  return _localized_vectors[&localized_vector].first;
}
const DenseVector< Number > & libMesh::DiffContext::get_localized_vector ( const NumericVector< Number > &  localized_vector) const [inherited]

const accessible version of get_localized_vector function

Definition at line 150 of file diff_context.C.

References libMesh::DiffContext::_localized_vectors, and libMesh::libmesh_assert().

{
  std::map<const NumericVector<Number>*, std::pair<DenseVector<Number>, std::vector<DenseSubVector<Number>*> > >::const_iterator
    localized_vectors_it = _localized_vectors.find(&localized_vector);
  libmesh_assert(localized_vectors_it != _localized_vectors.end());
  return localized_vectors_it->second.first;
}
const System* libMesh::FEMContext::get_mesh_system ( ) const [inline]

Accessor for moving mesh System

Definition at line 649 of file fem_context.h.

References _mesh_sys.

  { return this->_mesh_sys; }

Accessor for moving mesh System

Definition at line 655 of file fem_context.h.

References _mesh_sys.

  { return this->_mesh_sys; }
unsigned int libMesh::FEMContext::get_mesh_x_var ( ) const [inline]

Accessor for x-variable of moving mesh System

Definition at line 661 of file fem_context.h.

References _mesh_x_var.

Referenced by _do_elem_position_set(), and elem_position_get().

  { return _mesh_x_var; }
unsigned int libMesh::FEMContext::get_mesh_y_var ( ) const [inline]

Accessor for y-variable of moving mesh System

Definition at line 675 of file fem_context.h.

References _mesh_y_var.

Referenced by _do_elem_position_set(), and elem_position_get().

  { return _mesh_y_var; }
unsigned int libMesh::FEMContext::get_mesh_z_var ( ) const [inline]

Accessor for z-variable of moving mesh System

Definition at line 689 of file fem_context.h.

References _mesh_z_var.

Referenced by _do_elem_position_set(), and elem_position_get().

  { return _mesh_z_var; }
const std::vector<DenseVector<Number> >& libMesh::DiffContext::get_qoi_derivatives ( ) const [inline, inherited]

Const accessor for QoI derivatives.

Definition at line 299 of file diff_context.h.

References libMesh::DiffContext::_elem_qoi_derivative.

Referenced by pre_fe_reinit().

  { return _elem_qoi_derivative; }
std::vector<DenseVector<Number> >& libMesh::DiffContext::get_qoi_derivatives ( ) [inline, inherited]

Non-const accessor for QoI derivatives.

Definition at line 305 of file diff_context.h.

References libMesh::DiffContext::_elem_qoi_derivative.

  { return _elem_qoi_derivative; }
const DenseSubVector<Number>& libMesh::DiffContext::get_qoi_derivatives ( unsigned int  qoi,
unsigned int  var 
) const [inline, inherited]

Const accessor for QoI derivative of a particular qoi and variable corresponding to the index arguments.

Definition at line 312 of file diff_context.h.

References libMesh::DiffContext::_elem_qoi_subderivatives, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_qoi_subderivatives.size(), qoi);
    libmesh_assert_greater(_elem_qoi_subderivatives[qoi].size(), var);
    libmesh_assert(_elem_qoi_subderivatives[qoi][var]);
    return *(_elem_qoi_subderivatives[qoi][var]);
  }
DenseSubVector<Number>& libMesh::DiffContext::get_qoi_derivatives ( unsigned int  qoi,
unsigned int  var 
) [inline, inherited]

Non-const accessor for QoI derivative of a particular qoi and variable corresponding to the index arguments.

Definition at line 324 of file diff_context.h.

References libMesh::DiffContext::_elem_qoi_subderivatives, and libMesh::libmesh_assert().

  {
    libmesh_assert_greater(_elem_qoi_subderivatives.size(), qoi);
    libmesh_assert_greater(_elem_qoi_subderivatives[qoi].size(), var);
    libmesh_assert(_elem_qoi_subderivatives[qoi][var]);
    return *(_elem_qoi_subderivatives[qoi][var]);
  }
const std::vector<Number>& libMesh::DiffContext::get_qois ( ) const [inline, inherited]

Const accessor for QoI vector.

Definition at line 287 of file diff_context.h.

References libMesh::DiffContext::_elem_qoi.

  { return _elem_qoi; }
std::vector<Number>& libMesh::DiffContext::get_qois ( ) [inline, inherited]

Non-const accessor for QoI vector.

Definition at line 293 of file diff_context.h.

References libMesh::DiffContext::_elem_qoi.

  { return _elem_qoi; }
unsigned char libMesh::FEMContext::get_side ( ) const [inline]

Accessor for current side of Elem object

Definition at line 723 of file fem_context.h.

References side.

Referenced by side_fe_reinit().

  { return side; }
template<typename OutputShape >
void libMesh::FEMContext::get_side_fe ( unsigned int  var,
FEGenericBase< OutputShape > *&  fe 
) const [inline]

Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in the mesh. If you have lower dimensional elements in the mesh and need to query for those FE objects, use the alternative get_side_fe method.

Definition at line 261 of file fem_context.h.

References get_dim().

  { this->get_side_fe<OutputShape>(var,fe,this->get_dim()); }
FEBase* libMesh::FEMContext::get_side_fe ( unsigned int  var) const [inline]

Accessor for side finite element object for scalar-valued variable var for the largest dimension in the mesh. If you have lower dimensional elements in the mesh and need to query for those FE objects, use the alternative get_side_fe method.

Definition at line 270 of file fem_context.h.

References get_dim(), and get_side_fe().

Referenced by get_side_fe().

  { return this->get_side_fe(var,this->get_dim()); }
template<typename OutputShape >
void libMesh::FEMContext::get_side_fe ( unsigned int  var,
FEGenericBase< OutputShape > *&  fe,
unsigned char  dim 
) const [inline]

Accessor for edge/face (2D/3D) finite element object for variable var for dimension dim.

Definition at line 956 of file fem_context.h.

References _side_fe_var, and libMesh::libmesh_assert().

{
  libmesh_assert( !_side_fe_var[dim].empty() );
  libmesh_assert_less ( var, (_side_fe_var[dim].size() ) );
  fe = cast_ptr<FEGenericBase<OutputShape>*>( (_side_fe_var[dim][var] ) );
}
FEBase * libMesh::FEMContext::get_side_fe ( unsigned int  var,
unsigned char  dim 
) const [inline]

Accessor for side finite element object for scalar-valued variable var for dimension dim.

Definition at line 965 of file fem_context.h.

References _side_fe_var, and libMesh::libmesh_assert().

{
  libmesh_assert( !_side_fe_var[dim].empty() );
  libmesh_assert_less ( var, (_side_fe_var[dim].size() ) );
  return cast_ptr<FEBase*>( (_side_fe_var[dim][var] ) );
}
const QBase& libMesh::FEMContext::get_side_qrule ( ) const [inline]

Accessor for element side quadrature rule for the dimension of the current _elem.

Definition at line 612 of file fem_context.h.

References get_elem_dim(), and get_side_qrule().

Referenced by get_side_qrule().

  { return this->get_side_qrule(this->get_elem_dim()); }
const QBase& libMesh::FEMContext::get_side_qrule ( unsigned char  dim) const [inline]

Accessor for element side quadrature rule.

Definition at line 625 of file fem_context.h.

References _side_qrule, and libMesh::libmesh_assert().

  { libmesh_assert(_side_qrule[dim]);
    return *(this->_side_qrule[dim]); }
const System& libMesh::DiffContext::get_system ( ) const [inline, inherited]

Accessor for associated system.

Definition at line 107 of file diff_context.h.

References libMesh::DiffContext::_system.

Referenced by libMesh::DGFEMContext::neighbor_side_fe_reinit().

  { return _system; }
Real libMesh::DiffContext::get_system_time ( ) const [inline, inherited]

Accessor for the time variable stored in the system class.

Definition at line 367 of file diff_context.h.

References libMesh::DiffContext::system_time.

Referenced by _update_time_from_system().

  { return system_time; }
Real libMesh::DiffContext::get_time ( ) const [inline, inherited]

Accessor for the time for which the current nonlinear_solution is defined.

Definition at line 373 of file diff_context.h.

References libMesh::DiffContext::time.

  { return time; }
bool libMesh::FEMContext::has_elem ( ) const [inline]

Test for current Elem object

Definition at line 703 of file fem_context.h.

References _elem.

Referenced by build_new_fe(), elem_fe_reinit(), and pre_fe_reinit().

  { return (this->_elem != NULL); }

Reports if the boundary id is found on the current side

Definition at line 174 of file fem_context.C.

References _boundary_info, get_elem(), libMesh::BoundaryInfo::has_boundary_id(), and side.

{
  return _boundary_info.has_boundary_id(&(this->get_elem()), side, id);
}
template<typename OutputType >
template void libMesh::FEMContext::interior_curl< Gradient > ( unsigned int  var,
unsigned int  qp,
OutputType &  curl_u 
) const

Returns the curl of the solution variable var at the physical point p on the current element.

Definition at line 501 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, libMesh::FEGenericBase< OutputType >::get_curl_phi(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape> > &curl_phi = fe->get_curl_phi();

  // Accumulate solution curl
  curl_u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    curl_u.add_scaled(curl_phi[l][qp], coef(l));

  return;
}
template<typename OutputType >
template void libMesh::FEMContext::interior_div< Number > ( unsigned int  var,
unsigned int  qp,
OutputType &  div_u 
) const

Returns the divergence of the solution variable var at the physical point p on the current element.

Definition at line 534 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, libMesh::FEGenericBase< OutputType >::get_div_phi(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename
    TensorTools::IncrementRank
    <typename TensorTools::MakeReal<OutputType>::type>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence> > &div_phi = fe->get_div_phi();

  // Accumulate solution curl
  div_u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    div_u += div_phi[l][qp] * coef(l);

  return;
}
Gradient libMesh::FEMContext::interior_gradient ( unsigned int  var,
unsigned int  qp 
) const

Returns the gradient of the solution variable var at the quadrature point qp on the current element interior. This API currently present for backward compatibility.

Definition at line 374 of file fem_context.C.

{
  Gradient du;

  this->interior_gradient( var, qp, du );

  return du;
}
template<typename OutputType >
void libMesh::FEMContext::interior_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  du 
) const

Returns the gradient of the solution variable var at the quadrature point qp on the current element interior. This is the preferred API.

Definition at line 386 of file fem_context.C.

{
  this->some_interior_gradient <OutputType,&DiffContext::get_elem_solution>
    (var, qp, du);
}
template<typename OutputType >
template void libMesh::FEMContext::interior_gradients< Tensor > ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  interior_gradients_vector 
) const

Fills a vector with the gradient of the solution variable var at all the quadrature points in the current element interior. This is the preferred API.

Definition at line 397 of file fem_context.C.

References libMesh::FEGenericBase< OutputType >::get_dphi().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = fe->get_dphi();

  // Loop over all the q_points in this finite element
  for (unsigned int qp=0; qp != du_vals.size(); qp++)
    {
      OutputType &du = du_vals[qp];

      // Compute the gradient at this q_point
      du = 0;

      for (unsigned int l=0; l != n_dofs; l++)
        du.add_scaled(dphi[l][qp], coef(l));
    }

  return;
}
Tensor libMesh::FEMContext::interior_hessian ( unsigned int  var,
unsigned int  qp 
) const

Returns the hessian of the solution variable var at the quadrature point qp on the current element interior. This API currently present for backward compatibility.

Definition at line 436 of file fem_context.C.

{
  Tensor d2u;

  this->interior_hessian( var, qp, d2u );

  return d2u;
}
template<typename OutputType >
void libMesh::FEMContext::interior_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  d2u 
) const

Returns the hessian of the solution variable var at the quadrature point qp on the current element interior. This is the preferred API.

Definition at line 446 of file fem_context.C.

{
  this->some_interior_hessian <OutputType,&DiffContext::get_elem_solution>
    (var, qp, d2u);
}
template<typename OutputType >
template void libMesh::FEMContext::interior_hessians< Tensor > ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  d2u_vals 
) const

Fills a vector of hessians of the _system_vector at the all the quadrature points in the current element interior. This is the preferred API.

Definition at line 456 of file fem_context.C.

References libMesh::FEGenericBase< OutputType >::get_d2phi().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();

  // Loop over all the q_points in this finite element
  for (unsigned int qp=0; qp != d2u_vals.size(); qp++)
    {
      OutputType &d2u = d2u_vals[qp];

      // Compute the gradient at this q_point
      d2u = 0;

      for (unsigned int l=0; l != n_dofs; l++)
        d2u.add_scaled(d2phi[l][qp], coef(l));
    }

  return;
}
template<typename OutputType >
template void libMesh::FEMContext::interior_rate< Gradient > ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const

Returns the time derivative (rate) of the solution variable var at the quadrature point qp on the current element interior.

Definition at line 1364 of file fem_context.C.

{
  this->some_interior_value <OutputType,&DiffContext::get_elem_solution_rate>
    (var, qp, u);
}
Number libMesh::FEMContext::interior_value ( unsigned int  var,
unsigned int  qp 
) const

Returns the value of the solution variable var at the quadrature point qp on the current element interior. This API currently present for backward compatibility.

Definition at line 319 of file fem_context.C.

{
  Number u;

  this->interior_value( var, qp, u );

  return u;
}
template<typename OutputType >
void libMesh::FEMContext::interior_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const

Returns the value of the solution variable var at the quadrature point qp on the current element interior. This is the preferred API.

Definition at line 329 of file fem_context.C.

{
  this->some_interior_value <OutputType,&DiffContext::get_elem_solution>
    (var, qp, u);
}
template<typename OutputType >
template void libMesh::FEMContext::interior_values< Gradient > ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  interior_values_vector 
) const

Fills a vector of values of the _system_vector at the all the quadrature points in the current element interior.

Definition at line 338 of file fem_context.C.

References libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_localized_subvector(), and libMesh::FEGenericBase< OutputType >::get_phi().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);

  // Get the finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();

  // Loop over all the q_points on this element
  for (unsigned int qp=0; qp != u_vals.size(); qp++)
    {
      OutputType &u = u_vals[qp];

      // Compute the value at this q_point
      u = 0.;

      for (unsigned int l=0; l != n_dofs; l++)
        u += phi[l][qp] * coef(l);
    }

  return;
}
bool libMesh::DiffContext::is_adjoint ( ) const [inline, inherited]

Accessor for querying whether we need to do a primal or adjoint solve

Definition at line 413 of file diff_context.h.

References libMesh::DiffContext::_is_adjoint.

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

  { return _is_adjoint; }
bool& libMesh::DiffContext::is_adjoint ( ) [inline, inherited]

Accessor for setting whether we need to do a primal or adjoint solve

Definition at line 420 of file diff_context.h.

References libMesh::DiffContext::_is_adjoint.

  { return _is_adjoint; }
unsigned int libMesh::DiffContext::n_vars ( ) const [inline, inherited]

Number of variables in solution.

Definition at line 101 of file diff_context.h.

References libMesh::DiffContext::_dof_indices_var.

  { return cast_int<unsigned int>(_dof_indices_var.size()); }
void libMesh::FEMContext::nonlocal_reinit ( Real  theta) [virtual]

Gives derived classes the opportunity to reinitialize data needed for nonlocal calculations at a new point within a timestep

Reimplemented from libMesh::DiffContext.

Definition at line 1425 of file fem_context.C.

References _update_time_from_system(), and elem_fe_reinit().

{
  // Update the "time" variable of this context object
  this->_update_time_from_system(theta);

  // We can reuse the Elem FE safely here.
  elem_fe_reinit();
}
template<typename OutputType >
template void libMesh::FEMContext::point_curl< Gradient > ( unsigned int  var,
const Point p,
OutputType &  curl_u 
) const

Returns the curl of the solution variable var at the physical point p on the current element.

Definition at line 978 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape> >&  curl_phi = fe_new->get_curl_phi();

  curl_u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    curl_u.add_scaled(curl_phi[l][0], coef(l));

  return;
}
Gradient libMesh::FEMContext::point_gradient ( unsigned int  var,
const Point p 
) const

Returns the gradient of the solution variable var at the physical point p on the current element. This API currently present for backward compatibility.

Definition at line 875 of file fem_context.C.

{
  Gradient grad_u;

  this->point_gradient( var, p, grad_u );

  return grad_u;
}
template<typename OutputType >
void libMesh::FEMContext::point_gradient ( unsigned int  var,
const Point p,
OutputType &  grad_u 
) const

Returns the gradient of the solution variable var at the physical point p on the current element. This is the preferred API.

Definition at line 887 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >&  dphi = fe_new->get_dphi();

  grad_u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    grad_u.add_scaled(dphi[l][0], coef(l));

  return;
}
Tensor libMesh::FEMContext::point_hessian ( unsigned int  var,
const Point p 
) const

Returns the hessian of the solution variable var at the physical point p on the current element. This API currently present for backward compatibility.

Definition at line 926 of file fem_context.C.

{
  Tensor hess_u;

  this->point_hessian( var, p, hess_u );

  return hess_u;
}
template<typename OutputType >
void libMesh::FEMContext::point_hessian ( unsigned int  var,
const Point p,
OutputType &  hess_u 
) const

Returns the hessian of the solution variable var at the physical point p on the current element. This is the preferred API.

Definition at line 937 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >&  d2phi = fe_new->get_d2phi();

  hess_u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    hess_u.add_scaled(d2phi[l][0], coef(l));

  return;
}
template<typename OutputType >
void libMesh::FEMContext::point_rate ( unsigned int  var,
const Point p,
OutputType &  u 
) const

Returns the time derivative (rate) of the solution variable var at the physical point p on the current element.

Number libMesh::FEMContext::point_value ( unsigned int  var,
const Point p 
) const

Returns the value of the solution variable var at the physical point p on the current element. This API currently present for backward compatibility.

Definition at line 830 of file fem_context.C.

Referenced by libMesh::ParsedFEMFunction< Output >::component(), and libMesh::ParsedFEMFunction< Output >::operator()().

{
  Number u = 0.;

  this->point_value( var, p, u );

  return u;
}
template<typename OutputType >
void libMesh::FEMContext::point_value ( unsigned int  var,
const Point p,
OutputType &  u 
) const

Returns the value of the solution variable var at the physical point p on the current element. This is the preferred API.

Definition at line 840 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, build_new_fe(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Build a FE for calculating u(p)
  UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );

  // Get the values of the shape function derivatives
  const std::vector<std::vector<OutputShape> >&  phi = fe_new->get_phi();

  u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    u += phi[l][0] * coef(l);

  return;
}
void libMesh::FEMContext::pre_fe_reinit ( const System sys,
const Elem e 
) [virtual]

Reinitializes local data vectors/matrices on the current geometric element

Definition at line 1630 of file fem_context.C.

References libMesh::DiffContext::_elem_qoi_subderivatives, libMesh::DiffContext::_localized_vectors, libMesh::System::current_local_solution, libMesh::DofMap::dof_indices(), libMesh::NumericVector< T >::get(), libMesh::DiffContext::get_dof_indices(), libMesh::System::get_dof_map(), get_elem(), libMesh::DiffContext::get_elem_fixed_solution(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_residual(), libMesh::DiffContext::get_elem_solution(), libMesh::DiffContext::get_elem_solution_rate(), libMesh::DiffContext::get_qoi_derivatives(), libMesh::DenseVector< T >::get_values(), has_elem(), libMesh::System::n_vars(), libMesh::System::qoi, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), set_elem(), and libMesh::System::use_fixed_solution.

Referenced by libMesh::FEMSystem::assembly(), libMesh::FEMSystem::mesh_position_get(), libMesh::FEMSystem::mesh_position_set(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), and libMesh::System::project_vector().

{
  this->set_elem(e);

  // Initialize the per-element data for elem.
  if(this->has_elem())
    sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices());
  else
    // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
    sys.get_dof_map().dof_indices (NULL, this->get_dof_indices());

  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices().size());
  const std::size_t n_qoi = sys.qoi.size();

  // This also resizes elem_solution
  sys.current_local_solution->get(this->get_dof_indices(), this->get_elem_solution().get_values());

  if (sys.use_fixed_solution)
    this->get_elem_fixed_solution().resize(n_dofs);
  this->get_elem_solution_rate().resize(n_dofs);

  // These resize calls also zero out the residual and jacobian
  this->get_elem_residual().resize(n_dofs);
  this->get_elem_jacobian().resize(n_dofs, n_dofs);

  this->get_qoi_derivatives().resize(n_qoi);
  this->_elem_qoi_subderivatives.resize(n_qoi);
  for (std::size_t q=0; q != n_qoi; ++q)
    (this->get_qoi_derivatives())[q].resize(n_dofs);

  // Initialize the per-variable data for elem.
  {
    unsigned int sub_dofs = 0;
    for (unsigned int i=0; i != sys.n_vars(); ++i)
      {
        if(this->has_elem())
          sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);
        else
          // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
          sys.get_dof_map().dof_indices (NULL, this->get_dof_indices(i), i);


        const unsigned int n_dofs_var = cast_int<unsigned int>
          (this->get_dof_indices(i).size());

        this->get_elem_solution(i).reposition
          (sub_dofs, n_dofs_var);

        this->get_elem_solution_rate(i).reposition
          (sub_dofs, n_dofs_var);

        if (sys.use_fixed_solution)
          this->get_elem_fixed_solution(i).reposition
            (sub_dofs, n_dofs_var);

        this->get_elem_residual(i).reposition
          (sub_dofs, n_dofs_var);

        for (std::size_t q=0; q != n_qoi; ++q)
          this->get_qoi_derivatives(q,i).reposition
            (sub_dofs, n_dofs_var);

        for (unsigned int j=0; j != i; ++j)
          {
            const unsigned int n_dofs_var_j =
              cast_int<unsigned int>
              (this->get_dof_indices(j).size());

            this->get_elem_jacobian(i,j).reposition
              (sub_dofs, this->get_elem_residual(j).i_off(),
               n_dofs_var, n_dofs_var_j);
            this->get_elem_jacobian(j,i).reposition
              (this->get_elem_residual(j).i_off(), sub_dofs,
               n_dofs_var_j, n_dofs_var);
          }
        this->get_elem_jacobian(i,i).reposition
          (sub_dofs, sub_dofs,
           n_dofs_var,
           n_dofs_var);
        sub_dofs += n_dofs_var;
      }
    libmesh_assert_equal_to (sub_dofs, n_dofs);
  }

  // Now do the localization for the user requested vectors
  DiffContext::localized_vectors_iterator localized_vec_it = this->_localized_vectors.begin();
  const DiffContext::localized_vectors_iterator localized_vec_end = this->_localized_vectors.end();

  for(; localized_vec_it != localized_vec_end; ++localized_vec_it)
    {
      const NumericVector<Number>& current_localized_vector = *localized_vec_it->first;
      DenseVector<Number>& target_vector = localized_vec_it->second.first;

      current_localized_vector.get(this->get_dof_indices(), target_vector.get_values());

      // Initialize the per-variable data for elem.
      unsigned int sub_dofs = 0;
      for (unsigned int i=0; i != sys.n_vars(); ++i)
        {
          const unsigned int n_dofs_var = cast_int<unsigned int>
            (this->get_dof_indices(i).size());
          sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);

          localized_vec_it->second.second[i]->reposition
            (sub_dofs, n_dofs_var);

          sub_dofs += n_dofs_var;
        }
      libmesh_assert_equal_to (sub_dofs, n_dofs);
    }
}
void libMesh::DiffContext::set_deltat_pointer ( Real dt) [inherited]

Points the _deltat member of this class at a timestep value stored in the creating System, for example DiffSystem::deltat

Definition at line 113 of file diff_context.C.

References libMesh::DiffContext::_deltat.

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

{
  // We may actually want to be able to set this pointer to NULL, so
  // don't report an error for that.
  _deltat = dt;
}
void libMesh::FEMContext::set_elem ( const Elem e) [protected]

Helper function to promote accessor usage

Definition at line 1743 of file fem_context.C.

References _elem, _elem_dim, and libMesh::Elem::dim().

Referenced by pre_fe_reinit().

{
  this->_elem = e;

  // If e is NULL, we assume it's SCALAR and set _elem_dim to 0.
  this->_elem_dim = this->_elem ? this->_elem->dim() : 0;
}
virtual void libMesh::FEMContext::set_mesh_system ( System sys) [inline, virtual]

Tells the FEMContext 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.

This should be set automatically if the FEMPhysics requires it.

Definition at line 643 of file fem_context.h.

References _mesh_sys, and libMesh::sys.

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

  { this->_mesh_sys = sys; }
void libMesh::FEMContext::set_mesh_x_var ( unsigned int  x_var) [inline]

Accessor for x-variable of moving mesh System

This should be set automatically if the FEMPhysics requires it.

Definition at line 669 of file fem_context.h.

References _mesh_x_var.

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

  { _mesh_x_var = x_var; }
void libMesh::FEMContext::set_mesh_y_var ( unsigned int  y_var) [inline]

Accessor for y-variable of moving mesh System

This should be set automatically if the FEMPhysics requires it.

Definition at line 683 of file fem_context.h.

References _mesh_y_var.

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

  { _mesh_y_var = y_var; }
void libMesh::FEMContext::set_mesh_z_var ( unsigned int  z_var) [inline]

Accessor for z-variable of moving mesh System

This should be set automatically if the FEMPhysics requires it.

Definition at line 697 of file fem_context.h.

References _mesh_z_var.

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

  { _mesh_z_var = z_var; }
void libMesh::DiffContext::set_time ( Real  time_in) [inline, inherited]

Set the time for which the current nonlinear_solution is defined.

Definition at line 379 of file diff_context.h.

References libMesh::DiffContext::time.

Referenced by _update_time_from_system().

  { time = time_in; }

Lists the boundary ids found on the current side

Definition at line 180 of file fem_context.C.

References _boundary_info, libMesh::BoundaryInfo::boundary_ids(), get_elem(), and side.

{
  return _boundary_info.boundary_ids(&(this->get_elem()), side);
}

Reinitializes side FE objects on the current geometric element

Reimplemented in libMesh::DGFEMContext.

Definition at line 1458 of file fem_context.C.

References _side_fe, get_elem(), get_elem_dim(), get_side(), and libMesh::libmesh_assert().

Referenced by elem_side_reinit().

{
  // Initialize all the side FE objects on elem/side.
  // Logging of FE::reinit is done in the FE functions
  // We only reinit the FE objects for the current element
  // dimension
  const unsigned char dim = this->get_elem_dim();

  libmesh_assert( !_side_fe[dim].empty() );

  std::map<FEType, FEAbstract *>::iterator local_fe_end = _side_fe[dim].end();
  for (std::map<FEType, FEAbstract *>::iterator i = _side_fe[dim].begin();
       i != local_fe_end; ++i)
    {
      i->second->reinit(&(this->get_elem()), this->get_side());
    }
}
Gradient libMesh::FEMContext::side_gradient ( unsigned int  var,
unsigned int  qp 
) const

Returns the gradient of the solution variable var at the quadrature point qp on the current element side. This API currently present for backward compatibility.

Definition at line 649 of file fem_context.C.

{
  Gradient du;

  this->side_gradient( var, qp, du );

  return du;
}
template<typename OutputType >
void libMesh::FEMContext::side_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  du 
) const

Returns the gradient of the solution variable var at the quadrature point qp on the current element side. This is the preferred API.

Definition at line 660 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, libMesh::DiffContext::get_dof_indices(), libMesh::FEGenericBase< OutputType >::get_dphi(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector< typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();

  // Accumulate solution derivatives
  du = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    du.add_scaled(dphi[l][qp], coef(l));

  return;
}
template<typename OutputType >
template void libMesh::FEMContext::side_gradients< Tensor > ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  side_gradients_vector 
) const

Fills a vector with the gradient of the solution variable var at all the quadrature points on the current element side. This is the preferred API.

Definition at line 697 of file fem_context.C.

References libMesh::FEGenericBase< OutputType >::get_dphi().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();

  // Loop over all the q_points in this finite element
  for (unsigned int qp=0; qp != du_vals.size(); qp++)
    {
      OutputType &du = du_vals[qp];

      du = 0;

      // Compute the gradient at this q_point
      for (unsigned int l=0; l != n_dofs; l++)
        du.add_scaled(dphi[l][qp], coef(l));
    }

  return;
}
Tensor libMesh::FEMContext::side_hessian ( unsigned int  var,
unsigned int  qp 
) const

Returns the hessian of the solution variable var at the quadrature point qp on the current element side. This API currently present for backward compatibility.

Definition at line 736 of file fem_context.C.

{
  Tensor d2u;

  this->side_hessian( var, qp, d2u );

  return d2u;
}
template<typename OutputType >
void libMesh::FEMContext::side_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  d2u 
) const

Returns the hessian of the solution variable var at the quadrature point qp on the current element side. This is the preferred API.

Definition at line 746 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, libMesh::FEGenericBase< OutputType >::get_d2phi(), libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();

  // Accumulate solution second derivatives
  d2u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    d2u.add_scaled(d2phi[l][qp], coef(l));

  return;
}
template<typename OutputType >
template void libMesh::FEMContext::side_hessians< Tensor > ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  d2u_vals 
) const

Fills a vector of hessians of the _system_vector at the all the quadrature points on the current element side. This is the preferred API.

Definition at line 784 of file fem_context.C.

References libMesh::FEGenericBase< OutputType >::get_d2phi().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();

  // Loop over all the q_points in this finite element
  for (unsigned int qp=0; qp != d2u_vals.size(); qp++)
    {
      OutputType &d2u = d2u_vals[qp];

      // Compute the gradient at this q_point
      d2u = 0;

      for (unsigned int l=0; l != n_dofs; l++)
        d2u.add_scaled(d2phi[l][qp], coef(l));
    }

  return;
}
template<typename OutputType >
void libMesh::FEMContext::side_rate ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const

Returns the time derivative (rate) of the solution variable var at the quadrature point qp on the current element side.

Number libMesh::FEMContext::side_value ( unsigned int  var,
unsigned int  qp 
) const

Returns the value of the solution variable var at the quadrature point qp on the current element side. This API currently present for backward compatibility.

Definition at line 568 of file fem_context.C.

{
  Number u = 0.;

  this->side_value( var, qp, u );

  return u;
}
template<typename OutputType >
void libMesh::FEMContext::side_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const

Returns the value of the solution variable var at the quadrature point qp on the current element side. This is the preferred API.

Definition at line 579 of file fem_context.C.

References libMesh::DiffContext::_elem_subsolutions, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_solution(), libMesh::FEGenericBase< OutputType >::get_phi(), and libMesh::libmesh_assert().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  libmesh_assert_greater (this->_elem_subsolutions.size(), var);
  libmesh_assert(&(this->get_elem_solution(var)));
  const DenseSubVector<Number> &coef = this->get_elem_solution(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();

  // Accumulate solution value
  u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    u += phi[l][qp] * coef(l);

  return;
}
template<typename OutputType >
template void libMesh::FEMContext::side_values< Gradient > ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  side_values_vector 
) const

Fills a vector of values of the _system_vector at the all the quadrature points on the current element side.

Definition at line 613 of file fem_context.C.

References libMesh::FEGenericBase< OutputType >::get_phi().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);

  // Get the finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();

  // Loop over all the q_points on this element
  for (unsigned int qp=0; qp != u_vals.size(); qp++)
    {
      OutputType &u = u_vals[qp];

      // Compute the value at this q_point
      u = 0.;

      for (unsigned int l=0; l != n_dofs; l++)
        u += phi[l][qp] * coef(l);
    }

  return;
}
template<typename OutputType , FEMContext::diff_subsolution_getter subsolution_getter>
void libMesh::FEMContext::some_interior_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [protected]

Helper function to reduce some code duplication in the *interior_gradient methods.

Definition at line 219 of file fem_context.C.

References libMesh::DiffContext::get_dof_indices(), and libMesh::FEGenericBase< OutputType >::get_dphi().

{
  typedef typename TensorTools::MakeReal
    <typename TensorTools::DecrementRank<OutputType>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = fe->get_dphi();

  // Accumulate solution derivatives
  du = 0;

  for (unsigned int l=0; l != n_dofs; l++)
    du.add_scaled(dphi[l][qp], coef(l));

  return;
}
template<typename OutputType , FEMContext::diff_subsolution_getter subsolution_getter>
void libMesh::FEMContext::some_interior_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [protected]

Helper function to reduce some code duplication in the *interior_hessian methods.

Definition at line 254 of file fem_context.C.

References libMesh::FEGenericBase< OutputType >::get_d2phi(), and libMesh::DiffContext::get_dof_indices().

{
  typedef typename TensorTools::MakeReal<
    typename TensorTools::DecrementRank<
    typename TensorTools::DecrementRank<
    OutputType>::type>::type>::type
    OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();

  // Accumulate solution second derivatives
  d2u = 0.0;

  for (unsigned int l=0; l != n_dofs; l++)
    d2u.add_scaled(d2phi[l][qp], coef(l));

  return;
}
template<typename OutputType , FEMContext::diff_subsolution_getter subsolution_getter>
void libMesh::FEMContext::some_interior_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [protected]

Helper function to reduce some code duplication in the *interior_value methods.

Definition at line 189 of file fem_context.C.

References libMesh::DiffContext::get_dof_indices(), and libMesh::FEGenericBase< OutputType >::get_phi().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);

  // Get finite element object
  FEGenericBase<OutputShape>* fe = NULL;
  this->get_element_fe<OutputShape>( var, fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();

  // Accumulate solution value
  u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    u += phi[l][qp] * coef(l);
}
template<typename OutputType , FEMContext::diff_subsolution_getter subsolution_getter>
void libMesh::FEMContext::some_side_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [protected]

Helper function to reduce some code duplication in the *side_value methods.

Definition at line 291 of file fem_context.C.

References libMesh::DiffContext::get_dof_indices(), and libMesh::FEGenericBase< OutputType >::get_phi().

{
  typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;

  // Get local-to-global dof index lookup
  libmesh_assert_greater (this->get_dof_indices().size(), var);
  const unsigned int n_dofs = cast_int<unsigned int>
    (this->get_dof_indices(var).size());

  // Get current local coefficients
  const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);

  // Get finite element object
  FEGenericBase<OutputShape>* the_side_fe = NULL;
  this->get_side_fe<OutputShape>( var, the_side_fe );

  // Get shape function values at quadrature point
  const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();

  // Accumulate solution value
  u = 0.;

  for (unsigned int l=0; l != n_dofs; l++)
    u += phi[l][qp] * coef(l);

  return;
}

Member Data Documentation

Saved reference to BoundaryInfo on the mesh for this System. Used to answer boundary id requests.

Definition at line 863 of file fem_context.h.

Referenced by has_side_boundary_id(), and side_boundary_ids().

unsigned char libMesh::FEMContext::_dim [protected]

Cached dimension of largest dimension element in this mesh

Definition at line 873 of file fem_context.h.

Referenced by libMesh::DGFEMContext::DGFEMContext(), FEMContext(), and get_dim().

std::vector<dof_id_type> libMesh::DiffContext::_dof_indices [protected, inherited]

Global Degree of freedom index lists

Definition at line 566 of file diff_context.h.

Referenced by libMesh::DiffContext::get_dof_indices().

std::vector<std::vector<dof_id_type> > libMesh::DiffContext::_dof_indices_var [protected, inherited]

Definition at line 846 of file fem_context.h.

Referenced by edge_fe_reinit(), FEMContext(), and ~FEMContext().

std::vector<FEAbstract*> libMesh::FEMContext::_edge_fe_var [protected]

Definition at line 857 of file fem_context.h.

Referenced by FEMContext(), and get_edge_fe().

Quadrature rules for element edges. If the FEM context is told to prepare for edge integration on 3D elements, it will try to find a quadrature rule that correctly integrates all variables. Because edge rules only apply to 3D elements, we don't need to worry about multiple dimensions

Definition at line 903 of file fem_context.h.

Referenced by FEMContext(), get_edge_qrule(), and ~FEMContext().

const Elem* libMesh::FEMContext::_elem [protected]

Current element for element_* to examine

Definition at line 868 of file fem_context.h.

Referenced by get_elem(), has_elem(), and set_elem().

unsigned char libMesh::FEMContext::_elem_dim [protected]

Cached dimension of this->_elem.

Definition at line 878 of file fem_context.h.

Referenced by get_elem_dim(), and set_elem().

Element by element components of nonlinear_solution at a fixed point in a timestep, for optional use by e.g. stabilized methods

Definition at line 532 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_elem_fixed_solution().

Element jacobian: derivatives of elem_residual with respect to elem_solution

Definition at line 544 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_elem_jacobian().

std::vector<Number> libMesh::DiffContext::_elem_qoi [protected, inherited]

Element quantity of interest contributions

Definition at line 549 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_qois().

std::vector<DenseVector<Number> > libMesh::DiffContext::_elem_qoi_derivative [protected, inherited]

Element quantity of interest derivative contributions

Definition at line 554 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_qoi_derivatives().

Element residual vector

Definition at line 538 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_elem_residual().

Element by element components of nonlinear_solution as adjusted by a time_solver

Definition at line 517 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_elem_solution().

Element by element components of du/dt as adjusted by a time_solver

Definition at line 524 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), and libMesh::DiffContext::get_elem_solution_rate().

std::vector<DenseSubVector<Number> *> libMesh::DiffContext::_elem_subresiduals [protected, inherited]

Element residual subvectors and Jacobian submatrices

Definition at line 560 of file diff_context.h.

Referenced by libMesh::DiffContext::DiffContext(), libMesh::DiffContext::get_elem_residual(), and libMesh::DiffContext::~DiffContext().

std::vector<std::map<FEType, FEAbstract*> > libMesh::FEMContext::_element_fe [protected]

Finite element objects for each variable's interior, sides and edges. We store FE objects for each element dimension present in the mesh, except for edge_fe which only applies to 3D elements.

Definition at line 844 of file fem_context.h.

Referenced by elem_fe_reinit(), FEMContext(), and ~FEMContext().

std::vector<std::vector<FEAbstract*> > libMesh::FEMContext::_element_fe_var [protected]

Pointers to the same finite element objects, but indexed by variable number. We store FE objects for each element dimension present in the mesh, except for edge_fe_var which only applies for 3D elements.

Definition at line 855 of file fem_context.h.

Referenced by FEMContext(), and get_element_fe().

std::vector<QBase*> libMesh::FEMContext::_element_qrule [protected]

Quadrature rule for element interior. The FEM context will try to find a quadrature rule that correctly integrates all variables. We prepare quadrature rules for each element dimension in the mesh.

Definition at line 886 of file fem_context.h.

Referenced by FEMContext(), get_element_qrule(), and ~FEMContext().

std::map<const NumericVector<Number>*, std::pair<DenseVector<Number>, std::vector<DenseSubVector<Number>*> > > libMesh::DiffContext::_localized_vectors [protected, inherited]

Contains pointers to vectors the user has asked to be localized, keyed with pairs of element localized versions of that vector and per variable views

Definition at line 511 of file diff_context.h.

Referenced by libMesh::DiffContext::add_localized_vector(), libMesh::DiffContext::get_localized_subvector(), libMesh::DiffContext::get_localized_vector(), pre_fe_reinit(), and libMesh::DiffContext::~DiffContext().

Variables from which to acquire moving mesh information

Definition at line 769 of file fem_context.h.

Referenced by get_mesh_x_var(), and set_mesh_x_var().

Definition at line 769 of file fem_context.h.

Referenced by get_mesh_y_var(), and set_mesh_y_var().

Definition at line 769 of file fem_context.h.

Referenced by get_mesh_z_var(), and set_mesh_z_var().

std::vector<std::vector<FEAbstract*> > libMesh::FEMContext::_side_fe_var [protected]

Definition at line 856 of file fem_context.h.

Referenced by FEMContext(), and get_side_fe().

std::vector<QBase*> libMesh::FEMContext::_side_qrule [protected]

Quadrature rules for element sides The FEM context will try to find a quadrature rule that correctly integrates all variables. We prepare quadrature rules for each element dimension in the mesh.

Definition at line 894 of file fem_context.h.

Referenced by FEMContext(), get_side_qrule(), and ~FEMContext().

Current edge for edge_* to examine

Definition at line 779 of file fem_context.h.

Referenced by get_edge().

The derivative of elem_solution_rate with respect to the current nonlinear solution, for use by systems with non default mass_residual terms.

Definition at line 450 of file diff_context.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::EigenTimeSolver::element_residual(), libMesh::DiffContext::get_elem_solution_rate_derivative(), and libMesh::EigenTimeSolver::side_residual().

The derivative of elem_fixed_solution with respect to the nonlinear solution, for use by systems constructing jacobians with elem_fixed_solution based methods

Definition at line 457 of file diff_context.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), and libMesh::DiffContext::get_fixed_solution_derivative().

Current side for side_* to examine

Definition at line 774 of file fem_context.h.

Referenced by get_side(), has_side_boundary_id(), and side_boundary_ids().

This is the time stored in the System class at the time this context was created, i.e. the time at the beginning of the current timestep. This value gets set in the constructor and unlike DiffContext::time, is not tweaked mid-timestep by transient solvers: it remains equal to the value it was assigned at construction.

Definition at line 437 of file diff_context.h.

Referenced by libMesh::DiffContext::get_system_time().

For time-dependent problems, this is the time t for which the current nonlinear_solution is defined. FIXME - this needs to be tweaked mid-timestep by all transient solvers!

Definition at line 428 of file diff_context.h.

Referenced by libMesh::DiffContext::get_time(), and libMesh::DiffContext::set_time().


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