libMesh::FEMContext Class Reference

#include <fem_context.h>

Inheritance diagram for libMesh::FEMContext:

List of all members.

Public Types

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
template<typename OutputShape >
void get_side_fe (unsigned int var, FEGenericBase< OutputShape > *&fe) const
template<typename OutputShape >
void get_edge_fe (unsigned int var, FEGenericBase< OutputShape > *&fe) 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 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 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 pre_fe_reinit (const System &, const Elem *e)
void elem_fe_reinit ()
void side_fe_reinit ()
void edge_fe_reinit ()
const QBaseget_element_qrule () const
const QBaseget_side_qrule () 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)
const Elemget_elem () const
unsigned char get_side () const
unsigned char get_edge () const
unsigned char get_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
const DenseSubVector< Number > & get_elem_solution (unsigned int var) const
const DenseVector< Number > & get_elem_fixed_solution () const
const DenseSubVector< Number > & get_elem_fixed_solution (unsigned int var) const
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
const std::vector< dof_id_type > & get_dof_indices (unsigned int var) const
Real get_system_time () const
Real get_time () const
void set_time (Real time_in)
Real get_elem_solution_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

std::map< FEType, FEBase * > element_fe
std::map< FEType, FEBase * > side_fe
std::map< FEType, FEBase * > edge_fe
std::vector< FEBase * > element_fe_var
std::vector< FEBase * > side_fe_var
std::vector< FEBase * > edge_fe_var
QBaseelement_qrule
QBaseside_qrule
QBaseedge_qrule
System_mesh_sys
unsigned int _mesh_x_var
unsigned int _mesh_y_var
unsigned int _mesh_z_var
const Elemelem
unsigned char side
unsigned char edge
unsigned char dim
Real time
const Real system_time
DenseVector< Numberelem_solution
std::vector< DenseSubVector
< Number > * > 
elem_subsolutions
DenseVector< Numberelem_fixed_solution
std::vector< DenseSubVector
< Number > * > 
elem_fixed_subsolutions
Real elem_solution_derivative
Real fixed_solution_derivative
DenseVector< Numberelem_residual
DenseMatrix< Numberelem_jacobian
std::vector< Numberelem_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_typedof_indices
std::vector< std::vector
< dof_id_type > > 
dof_indices_var

Protected Member Functions

template<typename OutputShape >
AutoPtr< FEGenericBase
< OutputShape > > 
build_new_fe (const FEGenericBase< OutputShape > *fe, const Point &p) const

Protected Attributes

std::map< FEType, FEAbstract * > _element_fe
std::map< FEType, FEAbstract * > _side_fe
std::map< FEType, FEAbstract * > _edge_fe
std::vector< FEAbstract * > _element_fe_var
std::vector< FEAbstract * > _side_fe_var
std::vector< FEAbstract * > _edge_fe_var
BoundaryInfo_boundary_info
std::map< const NumericVector
< Number > *, std::pair
< DenseVector< Number >
, std::vector< DenseSubVector
< Number > * > > > 
localized_vectors

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 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 371 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 _edge_fe, _edge_fe_var, _element_fe, _element_fe_var, _side_fe, _side_fe_var, libMesh::FEAbstract::build(), libMesh::FEGenericBase< Real >::build(), libMesh::FEType::default_quadrature_rule(), dim, edge_fe, edge_fe_var, edge_qrule, element_fe, element_fe_var, element_qrule, libMesh::System::extra_quadrature_order, libMesh::FEInterface::field_type(), libMesh::System::n_vars(), libMesh::FEType::order, libMesh::AutoPtr< Tp >::release(), side_fe, side_fe_var, side_qrule, libMeshEnums::TYPE_SCALAR, and libMesh::System::variable_type().

00037   : DiffContext(sys),
00038     element_qrule(NULL), side_qrule(NULL),
00039     edge_qrule(NULL),
00040     _mesh_sys(NULL),
00041     _mesh_x_var(0),
00042     _mesh_y_var(0),
00043     _mesh_z_var(0),
00044     elem(NULL),
00045     side(0), edge(0), dim(sys.get_mesh().mesh_dimension()),
00046     _boundary_info(sys.get_mesh().boundary_info.get())
00047 {
00048   // We need to know which of our variables has the hardest
00049   // shape functions to numerically integrate.
00050 
00051   unsigned int nv = sys.n_vars();
00052 
00053   libmesh_assert (nv);
00054   FEType hardest_fe_type = sys.variable_type(0);
00055 
00056   for (unsigned int i=0; i != nv; ++i)
00057     {
00058       FEType fe_type = sys.variable_type(i);
00059 
00060       // FIXME - we don't yet handle mixed finite elements from
00061       // different families which require different quadrature rules
00062       // libmesh_assert_equal_to (fe_type.family, hardest_fe_type.family);
00063 
00064       if (fe_type.order > hardest_fe_type.order)
00065         hardest_fe_type = fe_type;
00066     }
00067 
00068   // Create an adequate quadrature rule
00069   element_qrule = hardest_fe_type.default_quadrature_rule
00070     (dim, sys.extra_quadrature_order).release();
00071   side_qrule = hardest_fe_type.default_quadrature_rule
00072     (dim-1, sys.extra_quadrature_order).release();
00073   if (dim == 3)
00074     edge_qrule = hardest_fe_type.default_quadrature_rule
00075       (1, sys.extra_quadrature_order).release();
00076 
00077   // Next, create finite element objects
00078   // Preserving backward compatibility here for now
00079   // Should move to the protected/FEAbstract interface
00080   element_fe_var.resize(nv);
00081   side_fe_var.resize(nv);
00082   if (dim == 3)
00083     edge_fe_var.resize(nv);
00084 
00085   _element_fe_var.resize(nv);
00086   _side_fe_var.resize(nv);
00087   if (dim == 3)
00088     _edge_fe_var.resize(nv);
00089 
00090   for (unsigned int i=0; i != nv; ++i)
00091     {
00092       FEType fe_type = sys.variable_type(i);
00093 
00094       // Preserving backward compatibility here for now
00095       // Should move to the protected/FEAbstract interface
00096       if( FEInterface::field_type( fe_type ) == TYPE_SCALAR )
00097         {
00098           if( element_fe[fe_type] == NULL )
00099             {
00100               element_fe[fe_type] = FEBase::build(dim, fe_type).release();
00101               element_fe[fe_type]->attach_quadrature_rule(element_qrule);
00102               side_fe[fe_type] = FEBase::build(dim, fe_type).release();
00103               side_fe[fe_type]->attach_quadrature_rule(side_qrule);
00104               
00105               if (dim == 3)
00106                 {
00107                   edge_fe[fe_type] = FEBase::build(dim, fe_type).release();
00108                   edge_fe[fe_type]->attach_quadrature_rule(edge_qrule);
00109                 }
00110             }
00111           element_fe_var[i] = element_fe[fe_type];
00112           side_fe_var[i] = side_fe[fe_type];
00113           if (dim == 3)
00114             edge_fe_var[i] = edge_fe[fe_type];
00115         }
00116 
00117       if ( _element_fe[fe_type] == NULL )
00118         {
00119           _element_fe[fe_type] = FEAbstract::build(dim, fe_type).release();
00120           _element_fe[fe_type]->attach_quadrature_rule(element_qrule);
00121           _side_fe[fe_type] = FEAbstract::build(dim, fe_type).release();
00122           _side_fe[fe_type]->attach_quadrature_rule(side_qrule);
00123 
00124           if (dim == 3)
00125             {
00126               _edge_fe[fe_type] = FEAbstract::build(dim, fe_type).release();
00127               _edge_fe[fe_type]->attach_quadrature_rule(edge_qrule);
00128             }
00129         }
00130       _element_fe_var[i] = _element_fe[fe_type];
00131       _side_fe_var[i] = _side_fe[fe_type];
00132       if (dim == 3)
00133         _edge_fe_var[i] = _edge_fe[fe_type];
00134       
00135     }
00136 }

libMesh::FEMContext::~FEMContext (  )  [virtual]

Destructor.

Definition at line 139 of file fem_context.C.

References _edge_fe, _element_fe, _side_fe, edge_fe, edge_qrule, element_fe, element_qrule, side_fe, and side_qrule.

00140 {
00141   // We don't want to store AutoPtrs in STL containers, but we don't
00142   // want to leak memory either
00143   // Preserving backward compatibility here for now
00144   // Should move to the protected/FEAbstract interface
00145   for (std::map<FEType, FEBase *>::iterator i = element_fe.begin();
00146        i != element_fe.end(); ++i)
00147     delete i->second;
00148   element_fe.clear();
00149 
00150   for (std::map<FEType, FEBase *>::iterator i = side_fe.begin();
00151        i != side_fe.end(); ++i)
00152     delete i->second;
00153   side_fe.clear();
00154 
00155   for (std::map<FEType, FEBase *>::iterator i = edge_fe.begin();
00156        i != edge_fe.end(); ++i)
00157     delete i->second;
00158   edge_fe.clear();
00159 
00160 
00161   for (std::map<FEType, FEAbstract *>::iterator i = _element_fe.begin();
00162        i != _element_fe.end(); ++i)
00163     delete i->second;
00164   _element_fe.clear();
00165   
00166   for (std::map<FEType, FEAbstract *>::iterator i = _side_fe.begin();
00167        i != _side_fe.end(); ++i)
00168     delete i->second;
00169   _side_fe.clear();
00170 
00171   for (std::map<FEType, FEAbstract *>::iterator i = _edge_fe.begin();
00172        i != _edge_fe.end(); ++i)
00173     delete i->second;
00174   _edge_fe.clear();
00175 
00176 
00177   delete element_qrule;
00178   element_qrule = NULL;
00179 
00180   delete side_qrule;
00181   side_qrule = NULL;
00182 
00183   if (edge_qrule)
00184     delete edge_qrule;
00185   side_qrule = NULL;
00186 }


Member Function Documentation

void libMesh::FEMContext::_do_elem_position_set ( Real  theta  )  [private]

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 1548 of file fem_context.C.

References _mesh_sys, _mesh_x_var, _mesh_y_var, _mesh_z_var, elem, libMesh::DiffContext::elem_subsolutions, element_fe_var, libMesh::invalid_uint, libMeshEnums::LAGRANGE, libMesh::libmesh_real(), libMesh::Elem::n_nodes(), n_nodes, and libMesh::n_threads().

Referenced by elem_position_set().

01549 {
01550   // This is too expensive to call unless we've been asked to move the mesh
01551   libmesh_assert (_mesh_sys);
01552 
01553   // This will probably break with threading when two contexts are
01554   // operating on elements which share a node
01555   libmesh_assert_equal_to (libMesh::n_threads(), 1);
01556 
01557   // If the coordinate data is in our own system, it's already
01558   // been set up for us, and we can ignore our input parameter theta
01559 //  if (_mesh_sys == this->number())
01560 //    {
01561       unsigned int n_nodes = elem->n_nodes();
01562       // For simplicity we demand that mesh coordinates be stored
01563       // in a format that allows a direct copy
01564       libmesh_assert(_mesh_x_var == libMesh::invalid_uint ||
01565                      (element_fe_var[_mesh_x_var]->get_fe_type().family
01566                       == LAGRANGE &&
01567                       elem_subsolutions[_mesh_x_var]->size() == n_nodes));
01568       libmesh_assert(_mesh_y_var == libMesh::invalid_uint ||
01569                      (element_fe_var[_mesh_y_var]->get_fe_type().family
01570                       == LAGRANGE &&
01571                       elem_subsolutions[_mesh_y_var]->size() == n_nodes));
01572       libmesh_assert(_mesh_z_var == libMesh::invalid_uint ||
01573                      (element_fe_var[_mesh_z_var]->get_fe_type().family
01574                       == LAGRANGE &&
01575                       elem_subsolutions[_mesh_z_var]->size() == n_nodes));
01576 
01577       // Set the new point coordinates
01578       if (_mesh_x_var != libMesh::invalid_uint)
01579         for (unsigned int i=0; i != n_nodes; ++i)
01580           const_cast<Elem*>(elem)->point(i)(0) =
01581             libmesh_real((*elem_subsolutions[_mesh_x_var])(i));
01582 
01583       if (_mesh_y_var != libMesh::invalid_uint)
01584         for (unsigned int i=0; i != n_nodes; ++i)
01585           const_cast<Elem*>(elem)->point(i)(1) =
01586             libmesh_real((*elem_subsolutions[_mesh_y_var])(i));
01587 
01588       if (_mesh_z_var != libMesh::invalid_uint)
01589         for (unsigned int i=0; i != n_nodes; ++i)
01590           const_cast<Elem*>(elem)->point(i)(2) =
01591             libmesh_real((*elem_subsolutions[_mesh_z_var])(i));
01592 //    }
01593   // FIXME - If the coordinate data is not in our own system, someone
01594   // had better get around to implementing that... - RHS
01595 //  else
01596 //    {
01597 //      libmesh_not_implemented();
01598 //    }
01599 }

void libMesh::FEMContext::_update_time_from_system ( Real  theta  )  [private]

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 1725 of file fem_context.C.

References libMesh::DiffContext::get_deltat_value(), libMesh::Real, libMesh::DiffContext::system_time, and libMesh::DiffContext::time.

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

01726 {
01727   // Update the "time" variable based on the value of theta.  For this
01728   // to work, we need to know the value of deltat, a pointer to which is now
01729   // stored by our parent DiffContext class.  Note: get_deltat_value() will
01730   // assert in debug mode if the requested pointer is NULL.
01731   const Real deltat = this->get_deltat_value();
01732 
01733   this->time = theta*(this->system_time + deltat) + (1.-theta)*this->system_time;
01734 }

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 125 of file diff_context.C.

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

00126 {
00127   // Make an empty pair keyed with a reference to this _localized_vector
00128   localized_vectors[&_localized_vector] = std::make_pair(DenseVector<Number>(), std::vector<DenseSubVector<Number>*>());
00129 
00130   unsigned int nv = _sys.n_vars();
00131 
00132   localized_vectors[&_localized_vector].second.reserve(nv);
00133   
00134   // Fill the DenseSubVector with nv copies of DenseVector
00135   for(unsigned int i=0; i != nv; ++i)
00136     localized_vectors[&_localized_vector].second.push_back(new DenseSubVector<Number>(localized_vectors[&_localized_vector].first));
00137 }

template<typename OutputShape >
AutoPtr< FEGenericBase< OutputShape > > libMesh::FEMContext::build_new_fe ( const FEGenericBase< OutputShape > *  fe,
const Point p 
) const [inline, protected]

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

Definition at line 1739 of file fem_context.C.

References dim, libMesh::Elem::dim(), elem, libMesh::FEAbstract::get_fe_type(), and libMesh::FEInterface::inverse_map().

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

01741 {
01742   FEType fe_type = fe->get_fe_type();
01743   AutoPtr<FEGenericBase<OutputShape> > fe_new(FEGenericBase<OutputShape>::build(elem->dim(), fe_type));
01744 
01745   // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h
01746   // Build a vector of point co-ordinates to send to reinit
01747   std::vector<Point> coor(1, FEInterface::inverse_map(dim, fe_type, elem, p));
01748 
01749   // Reinitialize the element and compute the shape function values at coor
01750   fe_new->reinit (elem, &coor);
01751 
01752   return fe_new;
01753 }

void libMesh::FEMContext::edge_fe_reinit (  ) 

Reinitializes edge FE objects on the current geometric element

Definition at line 1465 of file fem_context.C.

References _edge_fe, dim, edge, edge_fe, and elem.

Referenced by elem_edge_reinit().

01466 {
01467   libmesh_assert_equal_to (dim, 3);
01468 
01469   // Initialize all the interior FE objects on elem/edge.
01470   // Logging of FE::reinit is done in the FE functions
01471   std::map<FEType, FEBase *>::iterator fe_end = edge_fe.end();
01472   for (std::map<FEType, FEBase *>::iterator i = edge_fe.begin();
01473        i != fe_end; ++i)
01474     {
01475       i->second->edge_reinit(elem, edge);
01476     }
01477 
01478   std::map<FEType, FEAbstract *>::iterator local_fe_end = _edge_fe.end();
01479   for (std::map<FEType, FEAbstract *>::iterator i = _edge_fe.begin();
01480        i != local_fe_end; ++i)
01481     {
01482       i->second->edge_reinit(elem, edge);
01483     }
01484 }

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 1408 of file fem_context.C.

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

01409 {
01410   // Update the "time" variable of this context object
01411   this->_update_time_from_system(theta);
01412 
01413   // Handle a moving element if necessary
01414   if (_mesh_sys)
01415     {
01416       // FIXME - not threadsafe yet!
01417       elem_position_set(theta);
01418       edge_fe_reinit();
01419     }
01420 }

void libMesh::FEMContext::elem_fe_reinit (  ) 

Reinitializes interior FE objects on the current geometric element

Definition at line 1423 of file fem_context.C.

References _element_fe, elem, and element_fe.

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

01424 {
01425   // Initialize all the interior FE objects on elem.
01426   // Logging of FE::reinit is done in the FE functions
01427   std::map<FEType, FEBase *>::iterator fe_end = element_fe.end();
01428   for (std::map<FEType, FEBase *>::iterator i = element_fe.begin();
01429        i != fe_end; ++i)
01430     {
01431       i->second->reinit(elem);
01432     }
01433 
01434 
01435   std::map<FEType, FEAbstract *>::iterator local_fe_end = _element_fe.end();
01436   for (std::map<FEType, FEAbstract *>::iterator i = _element_fe.begin();
01437        i != local_fe_end; ++i)
01438     {
01439       i->second->reinit(elem);
01440     }
01441 }

void libMesh::FEMContext::elem_position_get (  ) 

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

Definition at line 1488 of file fem_context.C.

References _mesh_sys, _mesh_x_var, _mesh_y_var, _mesh_z_var, libMesh::Elem::default_order(), elem, libMesh::DiffContext::elem_subsolutions, element_fe_var, libMesh::invalid_uint, libMeshEnums::LAGRANGE, libMesh::Elem::n_nodes(), n_nodes, libMesh::n_threads(), and libMesh::Elem::point().

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

01489 {
01490   // This is too expensive to call unless we've been asked to move the mesh
01491   libmesh_assert (_mesh_sys);
01492 
01493   // This will probably break with threading when two contexts are
01494   // operating on elements which share a node
01495   libmesh_assert_equal_to (libMesh::n_threads(), 1);
01496 
01497   // If the coordinate data is in our own system, it's already
01498   // been set up for us
01499 //  if (_mesh_sys == this->number())
01500 //    {
01501       unsigned int n_nodes = elem->n_nodes();
01502       // For simplicity we demand that mesh coordinates be stored
01503       // in a format that allows a direct copy
01504       libmesh_assert(_mesh_x_var == libMesh::invalid_uint ||
01505                      (element_fe_var[_mesh_x_var]->get_fe_type().family
01506                       == LAGRANGE &&
01507                       element_fe_var[_mesh_x_var]->get_fe_type().order
01508                       == elem->default_order()));
01509       libmesh_assert(_mesh_y_var == libMesh::invalid_uint ||
01510                      (element_fe_var[_mesh_y_var]->get_fe_type().family
01511                       == LAGRANGE &&
01512                       element_fe_var[_mesh_y_var]->get_fe_type().order
01513                       == elem->default_order()));
01514       libmesh_assert(_mesh_z_var == libMesh::invalid_uint ||
01515                      (element_fe_var[_mesh_z_var]->get_fe_type().family
01516                       == LAGRANGE &&
01517                       element_fe_var[_mesh_z_var]->get_fe_type().order
01518                       == elem->default_order()));
01519 
01520       // Get degree of freedom coefficients from point coordinates
01521       if (_mesh_x_var != libMesh::invalid_uint)
01522         for (unsigned int i=0; i != n_nodes; ++i)
01523           (*elem_subsolutions[_mesh_x_var])(i) = elem->point(i)(0);
01524 
01525       if (_mesh_y_var != libMesh::invalid_uint)
01526         for (unsigned int i=0; i != n_nodes; ++i)
01527           (*elem_subsolutions[_mesh_y_var])(i) = elem->point(i)(1);
01528 
01529       if (_mesh_z_var != libMesh::invalid_uint)
01530         for (unsigned int i=0; i != n_nodes; ++i)
01531           (*elem_subsolutions[_mesh_z_var])(i) = elem->point(i)(2);
01532 //    }
01533   // FIXME - If the coordinate data is not in our own system, someone
01534   // had better get around to implementing that... - RHS
01535 //  else
01536 //    {
01537 //      libmesh_not_implemented();
01538 //    }
01539 }

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 745 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().

00746 {
00747   if (_mesh_sys)
00748     this->_do_elem_position_set(theta);
00749 }

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 1371 of file fem_context.C.

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

01372 {
01373   // Update the "time" variable of this context object
01374   this->_update_time_from_system(theta);
01375 
01376   // Handle a moving element if necessary.
01377   if (_mesh_sys)
01378     // We assume that the ``default'' state
01379     // of the mesh is its final, theta=1.0
01380     // position, so we don't bother with
01381     // mesh motion in that case.
01382     if (theta != 1.0)
01383       {
01384         // FIXME - ALE is not threadsafe yet!
01385         libmesh_assert_equal_to (libMesh::n_threads(), 1);
01386 
01387         elem_position_set(theta);
01388       }
01389   elem_fe_reinit();
01390 }

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 1393 of file fem_context.C.

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

01394 {
01395   // Update the "time" variable of this context object
01396   this->_update_time_from_system(theta);
01397 
01398   // Handle a moving element if necessary
01399   if (_mesh_sys)
01400     {
01401       // FIXME - not threadsafe yet!
01402       elem_position_set(theta);
01403       side_fe_reinit();
01404     }
01405 }

template<typename OutputType >
void libMesh::FEMContext::fixed_interior_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  grad_u 
) const [inline]

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

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 986 of file fem_context.C.

00987 {
00988   Gradient du;
00989 
00990   this->fixed_interior_gradient( var, qp, du );
00991 
00992   return du;
00993 }

template<typename OutputType >
void libMesh::FEMContext::fixed_interior_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  hess_u 
) const [inline]

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 1044 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_fixed_subsolutions, and libMesh::FEGenericBase< OutputType >::get_d2phi().

01046 {
01047   typedef typename TensorTools::MakeReal<
01048     typename TensorTools::DecrementRank<
01049       typename TensorTools::DecrementRank<
01050         OutputType>::type>::type>::type
01051     OutputShape;
01052 
01053   // Get local-to-global dof index lookup
01054   libmesh_assert_greater (dof_indices.size(), var);
01055   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01056     (dof_indices_var[var].size());
01057 
01058   // Get current local coefficients
01059   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01060   libmesh_assert(elem_fixed_subsolutions[var]);
01061   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01062 
01063   // Get finite element object
01064   FEGenericBase<OutputShape>* fe = NULL;
01065   this->get_element_fe<OutputShape>( var, fe );
01066 
01067   // Get shape function values at quadrature point
01068   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();
01069 
01070   // Accumulate solution second derivatives
01071   d2u = 0.0;
01072 
01073   for (unsigned int l=0; l != n_dofs; l++)
01074     d2u.add_scaled(d2phi[l][qp], coef(l));
01075 
01076   return;
01077 }

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 1033 of file fem_context.C.

01034 {
01035   Tensor d2u;
01036 
01037   this->fixed_interior_hessian( var, qp, d2u );
01038 
01039   return d2u;
01040 }

template<typename OutputType >
void libMesh::FEMContext::fixed_interior_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [inline]

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 953 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_fixed_subsolutions, and libMesh::FEGenericBase< OutputType >::get_phi().

00955 {
00956   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00957 
00958   // Get local-to-global dof index lookup
00959   libmesh_assert_greater (dof_indices.size(), var);
00960   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00961     (dof_indices_var[var].size());
00962 
00963   // Get current local coefficients
00964   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
00965   libmesh_assert(elem_fixed_subsolutions[var]);
00966   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
00967 
00968   // Get finite element object
00969   FEGenericBase<OutputShape>* fe = NULL;
00970   this->get_element_fe<OutputShape>( var, fe );
00971 
00972   // Get shape function values at quadrature point
00973   const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();
00974 
00975   // Accumulate solution value
00976   u = 0.0;
00977 
00978   for (unsigned int l=0; l != n_dofs; l++)
00979     u += phi[l][qp] * coef(l);
00980 
00981   return;
00982 }

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 941 of file fem_context.C.

00942 {
00943   Number u = 0.;
00944 
00945   this->fixed_interior_value( var, qp, u );
00946 
00947   return u;
00948 }

template<typename OutputType >
void libMesh::FEMContext::fixed_point_gradient ( unsigned int  var,
const Point p,
OutputType &  grad_u 
) const [inline]

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 1278 of file fem_context.C.

References build_new_fe(), libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, and libMesh::DiffContext::elem_fixed_subsolutions.

01280 {
01281   typedef typename TensorTools::MakeReal
01282     <typename TensorTools::DecrementRank<OutputType>::type>::type
01283       OutputShape;
01284 
01285   // Get local-to-global dof index lookup
01286   libmesh_assert_greater (dof_indices.size(), var);
01287   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01288     (dof_indices_var[var].size());
01289 
01290   // Get current local coefficients
01291   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01292   libmesh_assert(elem_fixed_subsolutions[var]);
01293   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01294 
01295   // Get finite element object
01296   FEGenericBase<OutputShape>* fe = NULL;
01297   this->get_element_fe<OutputShape>( var, fe );
01298 
01299   // Build a FE for calculating u(p)
01300   AutoPtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
01301 
01302   // Get the values of the shape function derivatives
01303   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >&  dphi = fe_new->get_dphi();
01304 
01305   grad_u = 0.0;
01306 
01307   for (unsigned int l=0; l != n_dofs; l++)
01308     grad_u.add_scaled(dphi[l][0], coef(l));
01309 
01310   return;
01311 }

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 1266 of file fem_context.C.

01267 {
01268   Gradient grad_u;
01269 
01270   this->fixed_point_gradient( var, p, grad_u );
01271 
01272   return grad_u;
01273 }

template<typename OutputType >
void libMesh::FEMContext::fixed_point_hessian ( unsigned int  var,
const Point p,
OutputType &  hess_u 
) const [inline]

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 1328 of file fem_context.C.

References build_new_fe(), libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, and libMesh::DiffContext::elem_fixed_subsolutions.

01330 {
01331   typedef typename TensorTools::MakeReal<
01332     typename TensorTools::DecrementRank<
01333       typename TensorTools::DecrementRank<
01334         OutputType>::type>::type>::type
01335     OutputShape;
01336 
01337   // Get local-to-global dof index lookup
01338   libmesh_assert_greater (dof_indices.size(), var);
01339   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01340     (dof_indices_var[var].size());
01341 
01342   // Get current local coefficients
01343   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01344   libmesh_assert(elem_fixed_subsolutions[var]);
01345   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01346 
01347   // Get finite element object
01348   FEGenericBase<OutputShape>* fe = NULL;
01349   this->get_element_fe<OutputShape>( var, fe );
01350   
01351   // Build a FE for calculating u(p)
01352   AutoPtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
01353 
01354   // Get the values of the shape function derivatives
01355   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >&  d2phi = fe_new->get_d2phi();
01356 
01357   hess_u = 0.0;
01358 
01359   for (unsigned int l=0; l != n_dofs; l++)
01360     hess_u.add_scaled(d2phi[l][0], coef(l));
01361 
01362   return;
01363 }

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 1316 of file fem_context.C.

01317 {
01318   Tensor hess_u;
01319 
01320   this->fixed_point_hessian( var, p, hess_u );
01321 
01322   return hess_u;
01323 }

template<typename OutputType >
void libMesh::FEMContext::fixed_point_value ( unsigned int  var,
const Point p,
OutputType &  u 
) const [inline]

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 1231 of file fem_context.C.

References build_new_fe(), libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, and libMesh::DiffContext::elem_fixed_subsolutions.

01233 {
01234   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
01235 
01236   // Get local-to-global dof index lookup
01237   libmesh_assert_greater (dof_indices.size(), var);
01238   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01239     (dof_indices_var[var].size());
01240 
01241   // Get current local coefficients
01242   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01243   libmesh_assert(elem_fixed_subsolutions[var]);
01244   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01245 
01246   // Get finite element object
01247   FEGenericBase<OutputShape>* fe = NULL;
01248   this->get_element_fe<OutputShape>( var, fe );
01249 
01250   // Build a FE for calculating u(p)
01251   AutoPtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
01252 
01253   // Get the values of the shape function derivatives
01254   const std::vector<std::vector<OutputShape> >&  phi = fe_new->get_phi();
01255   
01256   u = 0.;
01257 
01258   for (unsigned int l=0; l != n_dofs; l++)
01259     u += phi[l][0] * coef(l);
01260 
01261   return;
01262 }

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 1221 of file fem_context.C.

01222 {
01223   Number u = 0.;
01224 
01225   this->fixed_point_value( var, p, u );
01226 
01227   return u;
01228 }

template<typename OutputType >
void libMesh::FEMContext::fixed_side_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  grad_u 
) const [inline]

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 1137 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_fixed_subsolutions, and libMesh::FEGenericBase< OutputType >::get_dphi().

01139 {
01140   typedef typename TensorTools::MakeReal
01141     <typename TensorTools::DecrementRank<OutputType>::type>::type
01142       OutputShape;
01143 
01144   // Get local-to-global dof index lookup
01145   libmesh_assert_greater (dof_indices.size(), var);
01146   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01147     (dof_indices_var[var].size());
01148 
01149   // Get current local coefficients
01150   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01151   libmesh_assert(elem_fixed_subsolutions[var]);
01152   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01153 
01154   // Get finite element object
01155   FEGenericBase<OutputShape>* the_side_fe = NULL;
01156   this->get_side_fe<OutputShape>( var, the_side_fe );
01157 
01158   // Get shape function values at quadrature point
01159   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();
01160 
01161   // Accumulate solution derivatives
01162   du = 0.0;
01163 
01164   for (unsigned int l=0; l != n_dofs; l++)
01165     du.add_scaled(dphi[l][qp], coef(l));
01166 
01167   return;
01168 }

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 1126 of file fem_context.C.

01127 {
01128   Gradient du;
01129 
01130   this->fixed_side_gradient( var, qp, du );
01131 
01132   return du;
01133 }

template<typename OutputType >
void libMesh::FEMContext::fixed_side_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  hess_u 
) const [inline]

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 1183 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_fixed_subsolutions, and libMesh::FEGenericBase< OutputType >::get_d2phi().

01185 {
01186   typedef typename TensorTools::MakeReal<
01187     typename TensorTools::DecrementRank<
01188       typename TensorTools::DecrementRank<
01189         OutputType>::type>::type>::type
01190     OutputShape;
01191 
01192   // Get local-to-global dof index lookup
01193   libmesh_assert_greater (dof_indices.size(), var);
01194   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01195     (dof_indices_var[var].size());
01196 
01197   // Get current local coefficients
01198   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01199   libmesh_assert(elem_fixed_subsolutions[var]);
01200   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01201 
01202   // Get finite element object
01203   FEGenericBase<OutputShape>* the_side_fe = NULL;
01204   this->get_side_fe<OutputShape>( var, the_side_fe );
01205 
01206   // Get shape function values at quadrature point
01207   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();
01208 
01209   // Accumulate solution second derivatives
01210   d2u = 0.0;
01211 
01212   for (unsigned int l=0; l != n_dofs; l++)
01213     d2u.add_scaled(d2phi[l][qp], coef(l));
01214 
01215   return;
01216 }

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 1173 of file fem_context.C.

01174 {
01175   Tensor d2u;
01176 
01177   this->fixed_side_hessian( var, qp, d2u );
01178 
01179   return d2u;
01180 }

template<typename OutputType >
void libMesh::FEMContext::fixed_side_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [inline]

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 1093 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_fixed_subsolutions, and libMesh::FEGenericBase< OutputType >::get_phi().

01095 {
01096   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
01097 
01098   // Get local-to-global dof index lookup
01099   libmesh_assert_greater (dof_indices.size(), var);
01100   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01101     (dof_indices_var[var].size());
01102 
01103   // Get current local coefficients
01104   libmesh_assert_greater (elem_fixed_subsolutions.size(), var);
01105   libmesh_assert(elem_fixed_subsolutions[var]);
01106   DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
01107 
01108   // Get finite element object
01109   FEGenericBase<OutputShape>* the_side_fe = NULL;
01110   this->get_side_fe<OutputShape>( var, the_side_fe );
01111 
01112   // Get shape function values at quadrature point
01113   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
01114 
01115   // Accumulate solution value
01116   u = 0.0;
01117 
01118   for (unsigned int l=0; l != n_dofs; l++)
01119     u += phi[l][qp] * coef(l);
01120 
01121   return;
01122 }

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 1082 of file fem_context.C.

01083 {
01084   Number u = 0.;
01085 
01086   this->fixed_side_value( var, qp, u );
01087 
01088   return u;
01089 }

Real libMesh::DiffContext::get_deltat_value (  )  [inherited]

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

Definition at line 117 of file diff_context.C.

References libMesh::DiffContext::_deltat.

Referenced by _update_time_from_system().

00118 {
00119   libmesh_assert(_deltat);
00120 
00121   return *_deltat;
00122 }

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

Accessor for cached element dimension

Definition at line 607 of file fem_context.h.

References dim.

00608   { return dim; }

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 230 of file diff_context.h.

References libMesh::DiffContext::dof_indices_var.

00231   { return dof_indices_var[var]; }

const std::vector<dof_id_type>& libMesh::DiffContext::get_dof_indices (  )  const [inline, inherited]

Accessor for element dof indices

Definition at line 223 of file diff_context.h.

References libMesh::DiffContext::dof_indices.

00224   { return dof_indices; }

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

Accessor for current edge of Elem object

Definition at line 601 of file fem_context.h.

References edge.

00602   { 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 769 of file fem_context.h.

References _edge_fe_var.

00770 { 
00771   libmesh_assert_less ( var, _edge_fe_var.size() );
00772   fe = libmesh_cast_ptr<FEGenericBase<OutputShape>*>( _edge_fe_var[var] );
00773 }

const QBase* libMesh::FEMContext::get_edge_qrule (  )  const [inline]

Accessor for element edge quadrature rule.

Definition at line 518 of file fem_context.h.

References edge_qrule.

00519   { return this->edge_qrule; }

const Elem& libMesh::FEMContext::get_elem (  )  const [inline]

Accessor for current Elem object

Definition at line 589 of file fem_context.h.

References elem.

00590   { return *elem; }

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 127 of file diff_context.h.

References libMesh::DiffContext::elem_fixed_subsolutions.

00128   { return *(elem_fixed_subsolutions[var]); }

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

Accessor for element fixed solution.

Definition at line 120 of file diff_context.h.

References libMesh::DiffContext::elem_fixed_solution.

00121   { return elem_fixed_solution; }

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 179 of file diff_context.h.

References libMesh::DiffContext::elem_subjacobians.

00180   { return *(elem_subjacobians[var1][var2]); }

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

References libMesh::DiffContext::elem_subjacobians.

00173   { return *(elem_subjacobians[var1][var2]); }

DenseMatrix<Number>& libMesh::DiffContext::get_elem_jacobian (  )  [inline, inherited]

Non-const accessor for element Jacobian.

Definition at line 165 of file diff_context.h.

References libMesh::DiffContext::elem_jacobian.

00166   { return elem_jacobian; }

const DenseMatrix<Number>& libMesh::DiffContext::get_elem_jacobian (  )  const [inline, inherited]

Const accessor for element Jacobian.

Definition at line 159 of file diff_context.h.

References libMesh::DiffContext::elem_jacobian.

00160   { return elem_jacobian; }

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 153 of file diff_context.h.

References libMesh::DiffContext::elem_subresiduals.

00154   { return *(elem_subresiduals[var]); }

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 146 of file diff_context.h.

References libMesh::DiffContext::elem_subresiduals.

00147   { return *(elem_subresiduals[var]); }

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

Non-const accessor for element residual.

Definition at line 139 of file diff_context.h.

References libMesh::DiffContext::elem_residual.

00140   { return elem_residual; }

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

Const accessor for element residual.

Definition at line 133 of file diff_context.h.

References libMesh::DiffContext::elem_residual.

00134   { return elem_residual; }

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 114 of file diff_context.h.

References libMesh::DiffContext::elem_subsolutions.

00115   { return *(elem_subsolutions[var]); }

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

Accessor for element solution.

Definition at line 107 of file diff_context.h.

References libMesh::DiffContext::elem_solution.

00108   { return elem_solution; }

Real libMesh::DiffContext::get_elem_solution_derivative (  )  const [inline, inherited]

Definition at line 251 of file diff_context.h.

References libMesh::DiffContext::elem_solution_derivative.

00252   { return elem_solution_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.

Definition at line 753 of file fem_context.h.

References _element_fe_var.

00754 {
00755   libmesh_assert_less ( var, _element_fe_var.size() );
00756   fe = libmesh_cast_ptr<FEGenericBase<OutputShape>*>( _element_fe_var[var] );
00757 }

const QBase* libMesh::FEMContext::get_element_qrule (  )  const [inline]

Accessor for element interior quadrature rule.

Definition at line 506 of file fem_context.h.

References element_qrule.

00507   { return this->element_qrule; }

Real libMesh::DiffContext::get_fixed_solution_derivative (  )  const [inline, inherited]

Definition at line 254 of file diff_context.h.

References libMesh::DiffContext::fixed_solution_derivative.

00255   { return fixed_solution_derivative; }

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 161 of file diff_context.C.

References libMesh::DiffContext::localized_vectors.

00162 {
00163   std::map<const NumericVector<Number>*, std::pair<DenseVector<Number>, std::vector<DenseSubVector<Number>*> > >::const_iterator
00164     localized_vectors_it = localized_vectors.find(&_localized_vector);
00165   libmesh_assert(localized_vectors_it != localized_vectors.end());
00166   return *localized_vectors_it->second.second[_var];
00167 }

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 155 of file diff_context.C.

References libMesh::DiffContext::localized_vectors.

Referenced by interior_values().

00156 {
00157   return *localized_vectors[&_localized_vector].second[_var];
00158 }

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 146 of file diff_context.C.

References libMesh::DiffContext::localized_vectors.

00147 {
00148   std::map<const NumericVector<Number>*, std::pair<DenseVector<Number>, std::vector<DenseSubVector<Number>*> > >::const_iterator
00149     localized_vectors_it = localized_vectors.find(&_localized_vector);
00150   libmesh_assert(localized_vectors_it != localized_vectors.end());
00151   return localized_vectors_it->second.first;
00152 }

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 140 of file diff_context.C.

References libMesh::DiffContext::localized_vectors.

00141 {
00142   return localized_vectors[&_localized_vector].first;
00143 }

System* libMesh::FEMContext::get_mesh_system (  )  [inline]

Accessor for moving mesh System

Definition at line 541 of file fem_context.h.

References _mesh_sys.

00542   { return this->_mesh_sys; }

const System* libMesh::FEMContext::get_mesh_system (  )  const [inline]

Accessor for moving mesh System

Definition at line 535 of file fem_context.h.

References _mesh_sys.

00536   { 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 547 of file fem_context.h.

References _mesh_x_var.

00548   { 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 561 of file fem_context.h.

References _mesh_y_var.

00562   { 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 575 of file fem_context.h.

References _mesh_z_var.

00576   { return _mesh_z_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 217 of file diff_context.h.

References libMesh::DiffContext::elem_qoi_subderivatives.

00218   { return *(elem_qoi_subderivatives[qoi][var]); }

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 210 of file diff_context.h.

References libMesh::DiffContext::elem_qoi_subderivatives.

00211   { return *(elem_qoi_subderivatives[qoi][var]); }

std::vector<DenseVector<Number> >& libMesh::DiffContext::get_qoi_derivatives (  )  [inline, inherited]

Non-const accessor for QoI derivatives.

Definition at line 203 of file diff_context.h.

References libMesh::DiffContext::elem_qoi_derivative.

00204   { return elem_qoi_derivative; }

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

Const accessor for QoI derivatives.

Definition at line 197 of file diff_context.h.

References libMesh::DiffContext::elem_qoi_derivative.

00198   { return elem_qoi_derivative; }

std::vector<Number>& libMesh::DiffContext::get_qois (  )  [inline, inherited]

Non-const accessor for QoI vector.

Definition at line 191 of file diff_context.h.

References libMesh::DiffContext::elem_qoi.

00192   { return elem_qoi; }

const std::vector<Number>& libMesh::DiffContext::get_qois (  )  const [inline, inherited]

Const accessor for QoI vector.

Definition at line 185 of file diff_context.h.

References libMesh::DiffContext::elem_qoi.

00186   { return elem_qoi; }

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

Accessor for current side of Elem object

Definition at line 595 of file fem_context.h.

References side.

00596   { 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.

Definition at line 761 of file fem_context.h.

References _side_fe_var.

00762 {
00763   libmesh_assert_less ( var, _side_fe_var.size() );
00764   fe = libmesh_cast_ptr<FEGenericBase<OutputShape>*>( _side_fe_var[var] );
00765 }

const QBase* libMesh::FEMContext::get_side_qrule (  )  const [inline]

Accessor for element side quadrature rule.

Definition at line 512 of file fem_context.h.

References side_qrule.

00513   { return this->side_qrule; }

const System& libMesh::DiffContext::get_system (  )  const [inline, inherited]

Accessor for associated system.

Definition at line 101 of file diff_context.h.

References libMesh::DiffContext::_system.

00102   { return _system; }

Real libMesh::DiffContext::get_system_time (  )  const [inline, inherited]

Accessor for the time variable stored in the system class.

Definition at line 236 of file diff_context.h.

References libMesh::DiffContext::system_time.

00237   { 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 242 of file diff_context.h.

References libMesh::DiffContext::time.

00243   { return time; }

bool libMesh::FEMContext::has_side_boundary_id ( boundary_id_type  id  )  const

Reports if the boundary id is found on the current side

Definition at line 190 of file fem_context.C.

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

00191 {
00192   return _boundary_info->has_boundary_id(elem, side, id);
00193 }

template<typename OutputType >
void libMesh::FEMContext::interior_curl ( unsigned int  var,
unsigned int  qp,
OutputType &  curl_u 
) const [inline]

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

Definition at line 463 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_curl_phi().

00465 {
00466   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00467 
00468   // Get local-to-global dof index lookup
00469   libmesh_assert_greater (dof_indices.size(), var);
00470   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00471     (dof_indices_var[var].size());
00472 
00473   // Get current local coefficients
00474   libmesh_assert_greater (elem_subsolutions.size(), var);
00475   libmesh_assert(elem_subsolutions[var]);
00476   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00477 
00478   // Get finite element object
00479   FEGenericBase<OutputShape>* fe = NULL;
00480   this->get_element_fe<OutputShape>( var, fe );
00481 
00482   // Get shape function values at quadrature point
00483   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape> > &curl_phi = fe->get_curl_phi();
00484 
00485   // Accumulate solution curl
00486   curl_u = 0.;
00487 
00488   for (unsigned int l=0; l != n_dofs; l++)
00489     curl_u.add_scaled(curl_phi[l][qp], coef(l));
00490 
00491   return;
00492 }

template<typename OutputType >
void libMesh::FEMContext::interior_div ( unsigned int  var,
unsigned int  qp,
OutputType &  div_u 
) const [inline]

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

Definition at line 496 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_div_phi().

00498 {
00499   typedef typename 
00500     TensorTools::IncrementRank
00501       <typename TensorTools::MakeReal<OutputType>::type>::type OutputShape;
00502 
00503   // Get local-to-global dof index lookup
00504   libmesh_assert_greater (dof_indices.size(), var);
00505   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00506     (dof_indices_var[var].size());
00507 
00508   // Get current local coefficients
00509   libmesh_assert_greater (elem_subsolutions.size(), var);
00510   libmesh_assert(elem_subsolutions[var]);
00511   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00512 
00513   // Get finite element object
00514   FEGenericBase<OutputShape>* fe = NULL;
00515   this->get_element_fe<OutputShape>( var, fe );
00516 
00517   // Get shape function values at quadrature point
00518   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence> > &div_phi = fe->get_div_phi();
00519 
00520   // Accumulate solution curl
00521   div_u = 0.;
00522 
00523   for (unsigned int l=0; l != n_dofs; l++)
00524     div_u += div_phi[l][qp] * coef(l);
00525 
00526   return;
00527 }

template<typename OutputType >
void libMesh::FEMContext::interior_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  du 
) const [inline]

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 294 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_dphi().

00296 {
00297   typedef typename TensorTools::MakeReal
00298     <typename TensorTools::DecrementRank<OutputType>::type>::type
00299       OutputShape;
00300 
00301   // Get local-to-global dof index lookup
00302   libmesh_assert_greater (dof_indices.size(), var);
00303   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00304     (dof_indices_var[var].size());
00305 
00306   // Get current local coefficients
00307   libmesh_assert_greater (elem_subsolutions.size(), var);
00308   libmesh_assert(elem_subsolutions[var]);
00309   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00310 
00311   // Get finite element object
00312   FEGenericBase<OutputShape>* fe = NULL;
00313   this->get_element_fe<OutputShape>( var, fe );
00314 
00315   // Get shape function values at quadrature point
00316   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = fe->get_dphi();
00317 
00318   // Accumulate solution derivatives
00319   du = 0;
00320 
00321   for (unsigned int l=0; l != n_dofs; l++)
00322     du.add_scaled(dphi[l][qp], coef(l));
00323 
00324   return;
00325 }

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 282 of file fem_context.C.

00283 {
00284   Gradient du;
00285 
00286   this->interior_gradient( var, qp, du );
00287 
00288   return du;
00289 }

template<typename OutputType >
void libMesh::FEMContext::interior_gradients ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  interior_gradients_vector 
) const [inline]

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 331 of file fem_context.C.

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

00334 {
00335   typedef typename TensorTools::MakeReal
00336     <typename TensorTools::DecrementRank<OutputType>::type>::type
00337       OutputShape;
00338 
00339   // Get local-to-global dof index lookup
00340   libmesh_assert_greater (dof_indices.size(), var);
00341   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00342     (dof_indices_var[var].size());
00343 
00344   // Get current local coefficients  
00345   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00346 
00347   // Get finite element object
00348   FEGenericBase<OutputShape>* fe = NULL;
00349   this->get_element_fe<OutputShape>( var, fe );
00350 
00351   // Get shape function values at quadrature point
00352   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = fe->get_dphi();
00353   
00354   // Loop over all the q_points in this finite element
00355   for (unsigned int qp=0; qp != du_vals.size(); qp++)
00356     {      
00357       OutputType &du = du_vals[qp];
00358 
00359       // Compute the gradient at this q_point
00360       du = 0;
00361 
00362       for (unsigned int l=0; l != n_dofs; l++)
00363         du.add_scaled(dphi[l][qp], coef(l));
00364     }
00365 
00366   return;
00367 }

template<typename OutputType >
void libMesh::FEMContext::interior_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  d2u 
) const [inline]

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 380 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_d2phi().

00382 {
00383   typedef typename TensorTools::MakeReal<
00384     typename TensorTools::DecrementRank<
00385       typename TensorTools::DecrementRank<
00386         OutputType>::type>::type>::type
00387     OutputShape;
00388 
00389   // Get local-to-global dof index lookup
00390   libmesh_assert_greater (dof_indices.size(), var);
00391   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00392     (dof_indices_var[var].size());
00393 
00394   // Get current local coefficients
00395   libmesh_assert_greater (elem_subsolutions.size(), var);
00396   libmesh_assert(elem_subsolutions[var]);
00397   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00398 
00399   // Get finite element object
00400   FEGenericBase<OutputShape>* fe = NULL;
00401   this->get_element_fe<OutputShape>( var, fe );
00402 
00403   // Get shape function values at quadrature point
00404   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();
00405 
00406   // Accumulate solution second derivatives
00407   d2u = 0.0;
00408 
00409   for (unsigned int l=0; l != n_dofs; l++)
00410     d2u.add_scaled(d2phi[l][qp], coef(l));
00411 
00412   return;
00413 }

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 370 of file fem_context.C.

00371 {
00372   Tensor d2u;
00373 
00374   this->interior_hessian( var, qp, d2u );
00375 
00376   return d2u;
00377 }

template<typename OutputType >
void libMesh::FEMContext::interior_hessians ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  d2u_vals 
) const [inline]

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 418 of file fem_context.C.

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

00421 {
00422   typedef typename TensorTools::MakeReal<
00423     typename TensorTools::DecrementRank<
00424       typename TensorTools::DecrementRank<
00425         OutputType>::type>::type>::type
00426     OutputShape;
00427 
00428   // Get local-to-global dof index lookup
00429   libmesh_assert_greater (dof_indices.size(), var);
00430   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00431     (dof_indices_var[var].size());
00432 
00433   // Get current local coefficients  
00434   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00435 
00436   // Get finite element object
00437   FEGenericBase<OutputShape>* fe = NULL;
00438   this->get_element_fe<OutputShape>( var, fe );
00439 
00440   // Get shape function values at quadrature point
00441   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();
00442   
00443   // Loop over all the q_points in this finite element
00444   for (unsigned int qp=0; qp != d2u_vals.size(); qp++)
00445     {      
00446       OutputType &d2u = d2u_vals[qp];
00447 
00448       // Compute the gradient at this q_point
00449       d2u = 0;
00450 
00451       for (unsigned int l=0; l != n_dofs; l++)
00452         d2u.add_scaled(d2phi[l][qp], coef(l));
00453     }
00454 
00455   return;
00456 }

template<typename OutputType >
void libMesh::FEMContext::interior_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [inline]

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 213 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_phi().

00215 {
00216   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00217 
00218   // Get local-to-global dof index lookup
00219   libmesh_assert_greater (dof_indices.size(), var);
00220   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00221     (dof_indices_var[var].size());
00222 
00223   // Get current local coefficients
00224   libmesh_assert_greater (elem_subsolutions.size(), var);
00225   libmesh_assert(elem_subsolutions[var]);
00226   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00227 
00228   // Get finite element object
00229   FEGenericBase<OutputShape>* fe = NULL;
00230   this->get_element_fe<OutputShape>( var, fe );
00231 
00232   // Get shape function values at quadrature point
00233   const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();
00234 
00235   // Accumulate solution value
00236   u = 0.;
00237 
00238   for (unsigned int l=0; l != n_dofs; l++)
00239     u += phi[l][qp] * coef(l);
00240 
00241   return;
00242 }

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 203 of file fem_context.C.

00204 {
00205   Number u = 0.;
00206 
00207   this->interior_value( var, qp, u );
00208 
00209   return u;
00210 }

template<typename OutputType >
void libMesh::FEMContext::interior_values ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  interior_values_vector 
) const [inline]

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

Definition at line 246 of file fem_context.C.

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

00249 {
00250   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00251 
00252   // Get local-to-global dof index lookup
00253   libmesh_assert_greater (dof_indices.size(), var);
00254   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00255     (dof_indices_var[var].size());
00256   
00257   // Get current local coefficients  
00258   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00259       
00260   // Get the finite element object
00261   FEGenericBase<OutputShape>* fe = NULL;
00262   this->get_element_fe<OutputShape>( var, fe );
00263 
00264   // Get shape function values at quadrature point
00265   const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();
00266     
00267   // Loop over all the q_points on this element
00268   for (unsigned int qp=0; qp != u_vals.size(); qp++)
00269     {            
00270       OutputType &u = u_vals[qp];
00271       
00272       // Compute the value at this q_point
00273       u = 0.;
00274 
00275       for (unsigned int l=0; l != n_dofs; l++)
00276         u += phi[l][qp] * coef(l);
00277     }
00278   
00279   return;
00280 }

bool& libMesh::DiffContext::is_adjoint (  )  [inline, inherited]

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

Definition at line 268 of file diff_context.h.

References libMesh::DiffContext::_is_adjoint.

00269   { return _is_adjoint; }

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

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

Definition at line 261 of file diff_context.h.

References libMesh::DiffContext::_is_adjoint.

00262   { return _is_adjoint; }

unsigned int libMesh::DiffContext::n_vars (  )  const [inline, inherited]

Number of variables in solution.

Definition at line 95 of file diff_context.h.

References libMesh::DiffContext::dof_indices_var.

00096   { return libmesh_cast_int<unsigned int>(dof_indices_var.size()); }

template<typename OutputType >
void libMesh::FEMContext::point_gradient ( unsigned int  var,
const Point p,
OutputType &  grad_u 
) const [inline]

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 849 of file fem_context.C.

References build_new_fe(), libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, and libMesh::DiffContext::elem_subsolutions.

00851 {
00852   typedef typename TensorTools::MakeReal
00853     <typename TensorTools::DecrementRank<OutputType>::type>::type
00854       OutputShape;
00855 
00856   // Get local-to-global dof index lookup
00857   libmesh_assert_greater (dof_indices.size(), var);
00858   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00859     (dof_indices_var[var].size());
00860 
00861   // Get current local coefficients
00862   libmesh_assert_greater (elem_subsolutions.size(), var);
00863   libmesh_assert(elem_subsolutions[var]);
00864   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00865 
00866   // Get finite element object
00867   FEGenericBase<OutputShape>* fe = NULL;
00868   this->get_element_fe<OutputShape>( var, fe );
00869 
00870   // Build a FE for calculating u(p)
00871   AutoPtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00872 
00873   // Get the values of the shape function derivatives
00874   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >&  dphi = fe_new->get_dphi();
00875 
00876   grad_u = 0.0;
00877 
00878   for (unsigned int l=0; l != n_dofs; l++)
00879     grad_u.add_scaled(dphi[l][0], coef(l));
00880 
00881   return;
00882 }

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 837 of file fem_context.C.

00838 {
00839   Gradient grad_u;
00840 
00841   this->point_gradient( var, p, grad_u );
00842 
00843   return grad_u;
00844 }

template<typename OutputType >
void libMesh::FEMContext::point_hessian ( unsigned int  var,
const Point p,
OutputType &  hess_u 
) const [inline]

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 899 of file fem_context.C.

References build_new_fe(), libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, and libMesh::DiffContext::elem_subsolutions.

00901 {
00902   typedef typename TensorTools::MakeReal<
00903     typename TensorTools::DecrementRank<
00904       typename TensorTools::DecrementRank<
00905         OutputType>::type>::type>::type
00906     OutputShape;
00907 
00908   // Get local-to-global dof index lookup
00909   libmesh_assert_greater (dof_indices.size(), var);
00910   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00911     (dof_indices_var[var].size());
00912 
00913   // Get current local coefficients
00914   libmesh_assert_greater (elem_subsolutions.size(), var);
00915   libmesh_assert(elem_subsolutions[var]);
00916   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00917 
00918   // Get finite element object
00919   FEGenericBase<OutputShape>* fe = NULL;
00920   this->get_element_fe<OutputShape>( var, fe );
00921   
00922   // Build a FE for calculating u(p)
00923   AutoPtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00924 
00925   // Get the values of the shape function derivatives
00926   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >&  d2phi = fe_new->get_d2phi();
00927 
00928   hess_u = 0.0;
00929 
00930   for (unsigned int l=0; l != n_dofs; l++)
00931     hess_u.add_scaled(d2phi[l][0], coef(l));
00932 
00933   return;
00934 }

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 888 of file fem_context.C.

00889 {
00890   Tensor hess_u;
00891 
00892   this->point_hessian( var, p, hess_u );
00893 
00894   return hess_u;
00895 }

template<typename OutputType >
void libMesh::FEMContext::point_value ( unsigned int  var,
const Point p,
OutputType &  u 
) const [inline]

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 802 of file fem_context.C.

References build_new_fe(), libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, and libMesh::DiffContext::elem_subsolutions.

00804 {
00805   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00806 
00807   // Get local-to-global dof index lookup
00808   libmesh_assert_greater (dof_indices.size(), var);
00809   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00810     (dof_indices_var[var].size());
00811 
00812   // Get current local coefficients
00813   libmesh_assert_greater (elem_subsolutions.size(), var);
00814   libmesh_assert(elem_subsolutions[var]);
00815   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00816 
00817   // Get finite element object
00818   FEGenericBase<OutputShape>* fe = NULL;
00819   this->get_element_fe<OutputShape>( var, fe );
00820 
00821   // Build a FE for calculating u(p)
00822   AutoPtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00823 
00824   // Get the values of the shape function derivatives
00825   const std::vector<std::vector<OutputShape> >&  phi = fe_new->get_phi();
00826   
00827   u = 0.;
00828 
00829   for (unsigned int l=0; l != n_dofs; l++)
00830     u += phi[l][0] * coef(l);
00831 
00832   return;
00833 }

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 792 of file fem_context.C.

00793 {
00794   Number u = 0.;
00795 
00796   this->point_value( var, p, u );
00797 
00798   return u;
00799 }

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 1619 of file fem_context.C.

References libMesh::System::current_solution(), libMesh::DiffContext::dof_indices, libMesh::DofMap::dof_indices(), libMesh::DiffContext::dof_indices_var, elem, libMesh::DiffContext::elem_fixed_solution, libMesh::DiffContext::elem_fixed_subsolutions, libMesh::DiffContext::elem_jacobian, libMesh::DiffContext::elem_qoi_derivative, libMesh::DiffContext::elem_qoi_subderivatives, libMesh::DiffContext::elem_residual, libMesh::DiffContext::elem_solution, libMesh::DiffContext::elem_subjacobians, libMesh::DiffContext::elem_subresiduals, libMesh::DiffContext::elem_subsolutions, libMesh::System::get_dof_map(), libMesh::DiffContext::localized_vectors, libMesh::System::n_vars(), libMesh::System::qoi, libMesh::DenseMatrix< T >::resize(), libMesh::DenseVector< T >::resize(), and libMesh::System::use_fixed_solution.

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

01620 {
01621   elem = e;
01622 
01623   // Initialize the per-element data for elem.
01624   sys.get_dof_map().dof_indices (elem, dof_indices);
01625   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
01626     (dof_indices.size());
01627   const std::size_t n_qoi = sys.qoi.size();
01628 
01629   elem_solution.resize(n_dofs);
01630   if (sys.use_fixed_solution)
01631     elem_fixed_solution.resize(n_dofs);
01632 
01633   for (unsigned int i=0; i != n_dofs; ++i)
01634     elem_solution(i) = sys.current_solution(dof_indices[i]);
01635 
01636   // These resize calls also zero out the residual and jacobian
01637   elem_residual.resize(n_dofs);
01638   elem_jacobian.resize(n_dofs, n_dofs);
01639 
01640   elem_qoi_derivative.resize(n_qoi);
01641   elem_qoi_subderivatives.resize(n_qoi);
01642   for (std::size_t q=0; q != n_qoi; ++q)
01643     elem_qoi_derivative[q].resize(n_dofs);
01644 
01645   // Initialize the per-variable data for elem.
01646   {
01647     unsigned int sub_dofs = 0;
01648     for (unsigned int i=0; i != sys.n_vars(); ++i)
01649       {
01650         sys.get_dof_map().dof_indices (elem, dof_indices_var[i], i);
01651 
01652         const unsigned int n_dofs_var = libmesh_cast_int<unsigned int>
01653           (dof_indices_var[i].size());
01654 
01655         elem_subsolutions[i]->reposition
01656           (sub_dofs, n_dofs_var);
01657 
01658         if (sys.use_fixed_solution)
01659           elem_fixed_subsolutions[i]->reposition
01660             (sub_dofs, n_dofs_var);
01661 
01662         elem_subresiduals[i]->reposition
01663           (sub_dofs, n_dofs_var);
01664 
01665         for (std::size_t q=0; q != n_qoi; ++q)
01666           elem_qoi_subderivatives[q][i]->reposition
01667             (sub_dofs, n_dofs_var);
01668 
01669         for (unsigned int j=0; j != i; ++j)
01670           {
01671             const unsigned int n_dofs_var_j =
01672               libmesh_cast_int<unsigned int>
01673                 (dof_indices_var[j].size());
01674 
01675             elem_subjacobians[i][j]->reposition
01676               (sub_dofs, elem_subresiduals[j]->i_off(),
01677                n_dofs_var, n_dofs_var_j);
01678             elem_subjacobians[j][i]->reposition
01679               (elem_subresiduals[j]->i_off(), sub_dofs,
01680                n_dofs_var_j, n_dofs_var);
01681           }
01682         elem_subjacobians[i][i]->reposition
01683           (sub_dofs, sub_dofs,
01684            n_dofs_var,
01685            n_dofs_var);
01686         sub_dofs += n_dofs_var;
01687       }
01688     libmesh_assert_equal_to (sub_dofs, n_dofs);
01689   }
01690 
01691   // Now do the localization for the user requested vectors
01692   DiffContext::localized_vectors_iterator localized_vec_it = localized_vectors.begin();
01693   const DiffContext::localized_vectors_iterator localized_vec_end = localized_vectors.end();
01694   
01695   for(; localized_vec_it != localized_vec_end; ++localized_vec_it)
01696     {
01697       const NumericVector<Number> * current_localized_vector = localized_vec_it->first;
01698 
01699       DenseVector<Number>& target_vector = localized_vec_it->second.first;
01700       
01701       target_vector.resize(n_dofs);
01702           
01703       for (unsigned int i=0; i != n_dofs; ++i)
01704         target_vector(i) = (*current_localized_vector)(dof_indices[i]);
01705                           
01706       // Initialize the per-variable data for elem.
01707       unsigned int sub_dofs = 0;
01708       for (unsigned int i=0; i != sys.n_vars(); ++i)
01709         {
01710           const unsigned int n_dofs_var = libmesh_cast_int<unsigned int>
01711             (dof_indices_var[i].size());
01712           sys.get_dof_map().dof_indices (elem, dof_indices_var[i], i);
01713 
01714           localized_vec_it->second.second[i]->reposition
01715             (sub_dofs, n_dofs_var);
01716 
01717           sub_dofs += n_dofs_var;                     
01718         }
01719       libmesh_assert_equal_to (sub_dofs, n_dofs);         
01720     }
01721 }

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 109 of file diff_context.C.

References libMesh::DiffContext::_deltat.

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

00110 {
00111   // We may actually want to be able to set this pointer to NULL, so
00112   // don't report an error for that.
00113   _deltat = dt;
00114 }

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 529 of file fem_context.h.

References _mesh_sys.

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

00530   { 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 555 of file fem_context.h.

References _mesh_x_var.

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

00556   { _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 569 of file fem_context.h.

References _mesh_y_var.

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

00570   { _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 583 of file fem_context.h.

References _mesh_z_var.

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

00584   { _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 248 of file diff_context.h.

References libMesh::DiffContext::time.

00249   { time = time_in; }

std::vector< boundary_id_type > libMesh::FEMContext::side_boundary_ids (  )  const

Lists the boundary ids found on the current side

Definition at line 196 of file fem_context.C.

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

00197 {
00198   return _boundary_info->boundary_ids(elem, side);
00199 }

void libMesh::FEMContext::side_fe_reinit (  ) 

Reinitializes side FE objects on the current geometric element

Definition at line 1444 of file fem_context.C.

References _side_fe, elem, side, and side_fe.

Referenced by elem_side_reinit().

01445 {
01446   // Initialize all the side FE objects on elem/side.
01447   // Logging of FE::reinit is done in the FE functions
01448   std::map<FEType, FEBase *>::iterator fe_end = side_fe.end();
01449   for (std::map<FEType, FEBase *>::iterator i = side_fe.begin();
01450        i != fe_end; ++i)
01451     {
01452       i->second->reinit(elem, side);
01453     }
01454 
01455   std::map<FEType, FEAbstract *>::iterator local_fe_end = _side_fe.end();
01456   for (std::map<FEType, FEAbstract *>::iterator i = _side_fe.begin();
01457        i != local_fe_end; ++i)
01458     {
01459       i->second->reinit(elem, side);
01460     }
01461 }

template<typename OutputType >
void libMesh::FEMContext::side_gradient ( unsigned int  var,
unsigned int  qp,
OutputType &  du 
) const [inline]

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 622 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_dphi().

00624 {
00625   typedef typename TensorTools::MakeReal
00626     <typename TensorTools::DecrementRank<OutputType>::type>::type
00627       OutputShape;
00628 
00629   // Get local-to-global dof index lookup
00630   libmesh_assert_greater (dof_indices.size(), var);
00631   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00632     (dof_indices_var[var].size());
00633 
00634   // Get current local coefficients
00635   libmesh_assert_greater (elem_subsolutions.size(), var);
00636   libmesh_assert(elem_subsolutions[var]);
00637   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00638 
00639   // Get finite element object
00640   FEGenericBase<OutputShape>* the_side_fe = NULL;
00641   this->get_side_fe<OutputShape>( var, the_side_fe );
00642 
00643   // Get shape function values at quadrature point
00644   const std::vector<std::vector< typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();
00645 
00646   // Accumulate solution derivatives
00647   du = 0.;
00648 
00649   for (unsigned int l=0; l != n_dofs; l++)
00650     du.add_scaled(dphi[l][qp], coef(l));
00651 
00652   return;
00653 }

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 611 of file fem_context.C.

00612 {
00613   Gradient du;
00614 
00615   this->side_gradient( var, qp, du );
00616 
00617   return du;
00618 }

template<typename OutputType >
void libMesh::FEMContext::side_gradients ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  side_gradients_vector 
) const [inline]

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 659 of file fem_context.C.

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

00662 {
00663   typedef typename TensorTools::MakeReal
00664     <typename TensorTools::DecrementRank<OutputType>::type>::type
00665       OutputShape;
00666 
00667   // Get local-to-global dof index lookup
00668   libmesh_assert_greater (dof_indices.size(), var);
00669   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00670     (dof_indices_var[var].size());
00671 
00672   // Get current local coefficients  
00673   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00674 
00675   // Get finite element object
00676   FEGenericBase<OutputShape>* the_side_fe = NULL;
00677   this->get_side_fe<OutputShape>( var, the_side_fe );
00678 
00679   // Get shape function values at quadrature point
00680   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();
00681   
00682   // Loop over all the q_points in this finite element
00683   for (unsigned int qp=0; qp != du_vals.size(); qp++)
00684     {      
00685       OutputType &du = du_vals[qp];
00686 
00687       du = 0;
00688 
00689       // Compute the gradient at this q_point
00690       for (unsigned int l=0; l != n_dofs; l++)
00691         du.add_scaled(dphi[l][qp], coef(l));
00692     }
00693 
00694   return;
00695 }

template<typename OutputType >
void libMesh::FEMContext::side_hessian ( unsigned int  var,
unsigned int  qp,
OutputType &  d2u 
) const [inline]

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 708 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_d2phi().

00710 {
00711   typedef typename TensorTools::MakeReal<
00712     typename TensorTools::DecrementRank<
00713       typename TensorTools::DecrementRank<
00714         OutputType>::type>::type>::type
00715     OutputShape;
00716 
00717   // Get local-to-global dof index lookup
00718   libmesh_assert_greater (dof_indices.size(), var);
00719   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00720     (dof_indices_var[var].size());
00721 
00722   // Get current local coefficients
00723   libmesh_assert_greater (elem_subsolutions.size(), var);
00724   libmesh_assert(elem_subsolutions[var]);
00725   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00726 
00727   // Get finite element object
00728   FEGenericBase<OutputShape>* the_side_fe = NULL;
00729   this->get_side_fe<OutputShape>( var, the_side_fe );
00730 
00731   // Get shape function values at quadrature point
00732   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();
00733 
00734   // Accumulate solution second derivatives
00735   d2u = 0.0;
00736 
00737   for (unsigned int l=0; l != n_dofs; l++)
00738     d2u.add_scaled(d2phi[l][qp], coef(l));
00739 
00740   return;
00741 }

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 698 of file fem_context.C.

00699 {
00700   Tensor d2u;
00701 
00702   this->side_hessian( var, qp, d2u );
00703 
00704   return d2u;
00705 }

template<typename OutputType >
void libMesh::FEMContext::side_hessians ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  d2u_vals 
) const [inline]

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 746 of file fem_context.C.

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

00749 {
00750   typedef typename TensorTools::MakeReal<
00751     typename TensorTools::DecrementRank<
00752       typename TensorTools::DecrementRank<
00753         OutputType>::type>::type>::type
00754     OutputShape;
00755 
00756   // Get local-to-global dof index lookup
00757   libmesh_assert_greater (dof_indices.size(), var);
00758   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00759     (dof_indices_var[var].size());
00760 
00761   // Get current local coefficients  
00762   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00763 
00764   // Get finite element object
00765   FEGenericBase<OutputShape>* the_side_fe = NULL;
00766   this->get_side_fe<OutputShape>( var, the_side_fe );
00767 
00768   // Get shape function values at quadrature point
00769   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();
00770   
00771   // Loop over all the q_points in this finite element
00772   for (unsigned int qp=0; qp != d2u_vals.size(); qp++)
00773     {      
00774       OutputType &d2u = d2u_vals[qp];
00775 
00776       // Compute the gradient at this q_point
00777       d2u = 0;
00778 
00779       for (unsigned int l=0; l != n_dofs; l++)
00780         d2u.add_scaled(d2phi[l][qp], coef(l));
00781     }
00782 
00783   return;
00784 }

template<typename OutputType >
void libMesh::FEMContext::side_value ( unsigned int  var,
unsigned int  qp,
OutputType &  u 
) const [inline]

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 541 of file fem_context.C.

References libMesh::DiffContext::dof_indices, libMesh::DiffContext::dof_indices_var, libMesh::DiffContext::elem_subsolutions, and libMesh::FEGenericBase< OutputType >::get_phi().

00543 {
00544   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00545 
00546   // Get local-to-global dof index lookup
00547   libmesh_assert_greater (dof_indices.size(), var);
00548   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00549     (dof_indices_var[var].size());
00550 
00551   // Get current local coefficients
00552   libmesh_assert_greater (elem_subsolutions.size(), var);
00553   libmesh_assert(elem_subsolutions[var]);
00554   DenseSubVector<Number> &coef = *elem_subsolutions[var];
00555 
00556   // Get finite element object
00557   FEGenericBase<OutputShape>* the_side_fe = NULL;
00558   this->get_side_fe<OutputShape>( var, the_side_fe );
00559 
00560   // Get shape function values at quadrature point
00561   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
00562 
00563   // Accumulate solution value
00564   u = 0.;
00565 
00566   for (unsigned int l=0; l != n_dofs; l++)
00567     u += phi[l][qp] * coef(l);
00568 
00569   return;
00570 }

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 530 of file fem_context.C.

00531 {
00532   Number u = 0.;
00533 
00534   this->side_value( var, qp, u );
00535 
00536   return u;
00537 }

template<typename OutputType >
void libMesh::FEMContext::side_values ( unsigned int  var,
const NumericVector< Number > &  _system_vector,
std::vector< OutputType > &  side_values_vector 
) const [inline]

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

Definition at line 575 of file fem_context.C.

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

00578 {    
00579   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00580 
00581   // Get local-to-global dof index lookup
00582   libmesh_assert_greater (dof_indices.size(), var);
00583   const unsigned int n_dofs = libmesh_cast_int<unsigned int>
00584     (dof_indices_var[var].size());
00585   
00586   // Get current local coefficients  
00587   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00588       
00589   // Get the finite element object
00590   FEGenericBase<OutputShape>* the_side_fe = NULL;
00591   this->get_side_fe<OutputShape>( var, the_side_fe );
00592 
00593   // Get shape function values at quadrature point
00594   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
00595     
00596   // Loop over all the q_points on this element
00597   for (unsigned int qp=0; qp != u_vals.size(); qp++)
00598     {            
00599       OutputType &u = u_vals[qp];
00600 
00601       // Compute the value at this q_point
00602       u = 0.;
00603       
00604       for (unsigned int l=0; l != n_dofs; l++)
00605         u += phi[l][qp] * coef(l);
00606     }
00607   
00608   return;
00609 }


Member Data Documentation

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

Definition at line 718 of file fem_context.h.

Referenced by has_side_boundary_id(), and side_boundary_ids().

Definition at line 703 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 712 of file fem_context.h.

Referenced by FEMContext(), and get_edge_fe().

Finite element objects for each variable's interior, sides and edges.

Definition at line 701 of file fem_context.h.

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

Pointers to the same finite element objects, but indexed by variable number

Definition at line 710 of file fem_context.h.

Referenced by FEMContext(), and get_element_fe().

Variables from which to acquire moving mesh information

Definition at line 668 of file fem_context.h.

Referenced by _do_elem_position_set(), elem_position_get(), get_mesh_x_var(), and set_mesh_x_var().

Definition at line 702 of file fem_context.h.

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

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

Definition at line 711 of file fem_context.h.

Referenced by FEMContext(), and get_side_fe().

unsigned char libMesh::FEMContext::dim

Cached dimension of elements in this mesh

Definition at line 688 of file fem_context.h.

Referenced by build_new_fe(), edge_fe_reinit(), FEMContext(), and get_dim().

Current edge for edge_* to examine

Definition at line 683 of file fem_context.h.

Referenced by edge_fe_reinit(), and get_edge().

Definition at line 616 of file fem_context.h.

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

Definition at line 624 of file fem_context.h.

Referenced by FEMContext().

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

Definition at line 645 of file fem_context.h.

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

std::vector<Number> libMesh::DiffContext::elem_qoi [inherited]

Element quantity of interest contributions

Definition at line 330 of file diff_context.h.

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

Element quantity of interest derivative contributions

Definition at line 335 of file diff_context.h.

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

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

Definition at line 307 of file diff_context.h.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::DiffContext::get_elem_solution_derivative(), libMesh::EulerSolver::side_residual(), and libMesh::Euler2Solver::side_residual().

Element residual subvectors and Jacobian submatrices

Definition at line 341 of file diff_context.h.

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

Finite element objects for each variable's interior, sides and edges.

Definition at line 614 of file fem_context.h.

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

Pointers to the same finite element objects, but indexed by variable number

Definition at line 622 of file fem_context.h.

Referenced by _do_elem_position_set(), elem_position_get(), FEMContext(), and libMesh::FEMSystem::init_context().

Quadrature rule for element interior. The FEM context will try to find a quadrature rule that correctly integrates all variables

Definition at line 631 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 427 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().

Current side for side_* to examine

Definition at line 678 of file fem_context.h.

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

Definition at line 615 of file fem_context.h.

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

Definition at line 623 of file fem_context.h.

Referenced by FEMContext().

Quadrature rules for element sides The FEM context will try to find a quadrature rule that correctly integrates all variables

Definition at line 638 of file fem_context.h.

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

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 285 of file diff_context.h.

Referenced by _update_time_from_system(), and 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 276 of file diff_context.h.

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


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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:25 UTC

Hosted By:
SourceForge.net Logo