#include <fem_context.h>
FEMContext (const FEMSystem &sys)
virtual ~FEMContext ()
Number interior_value (unsigned int var, unsigned int qp)
Number side_value (unsigned int var, unsigned int qp)
Number point_value (unsigned int var, Point &p)
Gradient interior_gradient (unsigned int var, unsigned int qp)
Gradient side_gradient (unsigned int var, unsigned int qp)
Tensor interior_hessian (unsigned int var, unsigned int qp)
Tensor side_hessian (unsigned int var, unsigned int qp)
Number fixed_interior_value (unsigned int var, unsigned int qp)
Number fixed_side_value (unsigned int var, unsigned int qp)
Number fixed_point_value (unsigned int var, Point &p)
Gradient fixed_interior_gradient (unsigned int var, unsigned int qp)
Gradient fixed_side_gradient (unsigned int var, unsigned int qp)
Tensor fixed_interior_hessian (unsigned int var, unsigned int qp)
Tensor fixed_side_hessian (unsigned int var, unsigned int qp)
virtual void reinit (const FEMSystem &, Elem *)
virtual void elem_reinit (Real theta)
virtual void elem_side_reinit (Real theta)
void elem_fe_reinit ()
void elem_side_fe_reinit ()
void elem_position_set (Real theta)
void elem_position_get ()
std::map< FEType, FEBase * > element_fe
std::map< FEType, FEBase * > side_fe
std::vector< FEBase * > element_fe_var
std::vector< FEBase * > side_fe_var
QBase * element_qrule
QBase * side_qrule
unsigned int _mesh_sys
unsigned int _mesh_x_var
unsigned int _mesh_y_var
unsigned int _mesh_z_var
Elem * elem
unsigned char side
unsigned char dim
bool postprocess_sides
Real time
DenseVector< Number > elem_solution
std::vector< DenseSubVector< Number > * > elem_subsolutions
DenseVector< Number > elem_fixed_solution
std::vector< DenseSubVector< Number > * > elem_fixed_subsolutions
bool request_jacobian
Real elem_solution_derivative
Real fixed_solution_derivative
DenseVector< Number > elem_residual
DenseMatrix< Number > elem_jacobian
std::vector< DenseSubVector< Number > * > elem_subresiduals
std::vector< std::vector< DenseSubMatrix< Number > * > > elem_subjacobians
std::vector< unsigned int > dof_indices
std::vector< std::vector< unsigned int > > dof_indices_var
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:
Definition at line 58 of file fem_context.h.
Definition at line 41 of file fem_context.C.
References FEBase::build(), FEType::default_quadrature_rule(), dim, element_fe, element_fe_var, element_qrule, FEMSystem::extra_quadrature_order, System::n_vars(), FEType::order, AutoPtr< Tp >::release(), side_fe, side_fe_var, side_qrule, and System::variable_type().
: DiffContext(sys),
element_qrule(NULL), side_qrule(NULL),
_mesh_sys(libMesh::invalid_uint),
_mesh_x_var(libMesh::invalid_uint),
_mesh_y_var(libMesh::invalid_uint),
_mesh_z_var(libMesh::invalid_uint),
elem(NULL), side(0), dim(sys.get_mesh().mesh_dimension())
{
// We need to know which of our variables has the hardest
// shape functions to numerically integrate.
unsigned int n_vars = sys.n_vars();
libmesh_assert (n_vars);
FEType hardest_fe_type = sys.variable_type(0);
for (unsigned int i=0; i != n_vars; ++i)
{
FEType fe_type = sys.variable_type(i);
// FIXME - we don't yet handle mixed finite elements from
// different families which require different quadrature rules
// libmesh_assert (fe_type.family == hardest_fe_type.family);
if (fe_type.order > hardest_fe_type.order)
hardest_fe_type = fe_type;
}
// Create an adequate quadrature rule
element_qrule = hardest_fe_type.default_quadrature_rule
(dim, sys.extra_quadrature_order).release();
side_qrule = hardest_fe_type.default_quadrature_rule
(dim-1, sys.extra_quadrature_order).release();
// Next, create finite element objects
element_fe_var.resize(n_vars);
side_fe_var.resize(n_vars);
for (unsigned int i=0; i != n_vars; ++i)
{
FEType fe_type = sys.variable_type(i);
if (element_fe[fe_type] == NULL)
{
element_fe[fe_type] = FEBase::build(dim, fe_type).release();
element_fe[fe_type]->attach_quadrature_rule(element_qrule);
side_fe[fe_type] = FEBase::build(dim, fe_type).release();
side_fe[fe_type]->attach_quadrature_rule(side_qrule);
}
element_fe_var[i] = element_fe[fe_type];
side_fe_var[i] = side_fe[fe_type];
}
}
Definition at line 95 of file fem_context.C.
References element_fe, element_qrule, side_fe, and side_qrule.
{
// We don't want to store AutoPtrs in STL containers, but we don't
// want to leak memory either
for (std::map<FEType, FEBase *>::iterator i = element_fe.begin();
i != element_fe.end(); ++i)
delete i->second;
element_fe.clear();
for (std::map<FEType, FEBase *>::iterator i = side_fe.begin();
i != side_fe.end(); ++i)
delete i->second;
side_fe.clear();
delete element_qrule;
element_qrule = NULL;
delete side_qrule;
side_qrule = NULL;
}
Definition at line 514 of file fem_context.C.
References elem, and element_fe.
Referenced by elem_reinit(), and reinit().
{
// Initialize all the interior FE objects on elem.
// Logging of FE::reinit is done in the FE functions
std::map<FEType, FEBase *>::iterator fe_end = element_fe.end();
for (std::map<FEType, FEBase *>::iterator i = element_fe.begin();
i != fe_end; ++i)
{
i->second->reinit(elem);
}
}
Definition at line 541 of file fem_context.C.
References _mesh_sys, _mesh_x_var, _mesh_y_var, _mesh_z_var, Elem::default_order(), elem, DiffContext::elem_subsolutions, element_fe_var, libMesh::invalid_uint, libMeshEnums::LAGRANGE, Elem::n_nodes(), MeshTools::n_nodes(), and Elem::point().
Referenced by FEMSystem::mesh_position_get().
{
// This is too expensive to call unless we've been asked to move the mesh
libmesh_assert (_mesh_sys != libMesh::invalid_uint);
// This will probably break with threading when two contexts are
// operating on elements which share a node
libmesh_experimental();
// If the coordinate data is in our own system, it's already
// been set up for us
// if (_mesh_sys == this->number())
// {
unsigned int n_nodes = elem->n_nodes();
// For simplicity we demand that mesh coordinates be stored
// in a format that allows a direct copy
libmesh_assert(_mesh_x_var == libMesh::invalid_uint ||
(element_fe_var[_mesh_x_var]->get_fe_type().family
== LAGRANGE &&
element_fe_var[_mesh_x_var]->get_fe_type().order
== elem->default_order()));
libmesh_assert(_mesh_y_var == libMesh::invalid_uint ||
(element_fe_var[_mesh_y_var]->get_fe_type().family
== LAGRANGE &&
element_fe_var[_mesh_y_var]->get_fe_type().order
== elem->default_order()));
libmesh_assert(_mesh_z_var == libMesh::invalid_uint ||
(element_fe_var[_mesh_z_var]->get_fe_type().family
== LAGRANGE &&
element_fe_var[_mesh_z_var]->get_fe_type().order
== elem->default_order()));
// Get degree of freedom coefficients from point coordinates
if (_mesh_x_var != libMesh::invalid_uint)
for (unsigned int i=0; i != n_nodes; ++i)
(*elem_subsolutions[_mesh_x_var])(i) = elem->point(i)(0);
if (_mesh_y_var != libMesh::invalid_uint)
for (unsigned int i=0; i != n_nodes; ++i)
(*elem_subsolutions[_mesh_y_var])(i) = elem->point(i)(1);
if (_mesh_z_var != libMesh::invalid_uint)
for (unsigned int i=0; i != n_nodes; ++i)
(*elem_subsolutions[_mesh_z_var])(i) = elem->point(i)(2);
// }
// FIXME - If the coordinate data is not in our own system, someone
// had better get around to implementing that... - RHS
// else
// {
// libmesh_not_implemented();
// }
}
Definition at line 596 of file fem_context.C.
References _mesh_sys, _mesh_x_var, _mesh_y_var, _mesh_z_var, elem, DiffContext::elem_subsolutions, element_fe_var, libMesh::invalid_uint, libMeshEnums::LAGRANGE, libmesh_real(), Elem::n_nodes(), MeshTools::n_nodes(), and Elem::point().
Referenced by elem_reinit(), elem_side_reinit(), FEMSystem::mesh_position_set(), and reinit().
{
// This is too expensive to call unless we've been asked to move the mesh
libmesh_assert (_mesh_sys != libMesh::invalid_uint);
// This will probably break with threading when two contexts are
// operating on elements which share a node
libmesh_experimental();
// If the coordinate data is in our own system, it's already
// been set up for us, and we can ignore our input parameter theta
// if (_mesh_sys == this->number())
// {
unsigned int n_nodes = elem->n_nodes();
// For simplicity we demand that mesh coordinates be stored
// in a format that allows a direct copy
libmesh_assert(_mesh_x_var == libMesh::invalid_uint ||
(element_fe_var[_mesh_x_var]->get_fe_type().family
== LAGRANGE &&
elem_subsolutions[_mesh_x_var]->size() == n_nodes));
libmesh_assert(_mesh_y_var == libMesh::invalid_uint ||
(element_fe_var[_mesh_y_var]->get_fe_type().family
== LAGRANGE &&
elem_subsolutions[_mesh_y_var]->size() == n_nodes));
libmesh_assert(_mesh_z_var == libMesh::invalid_uint ||
(element_fe_var[_mesh_z_var]->get_fe_type().family
== LAGRANGE &&
elem_subsolutions[_mesh_z_var]->size() == n_nodes));
// Set the new point coordinates
if (_mesh_x_var != libMesh::invalid_uint)
for (unsigned int i=0; i != n_nodes; ++i)
elem->point(i)(0) =
libmesh_real((*elem_subsolutions[_mesh_x_var])(i));
if (_mesh_y_var != libMesh::invalid_uint)
for (unsigned int i=0; i != n_nodes; ++i)
elem->point(i)(1) =
libmesh_real((*elem_subsolutions[_mesh_y_var])(i));
if (_mesh_z_var != libMesh::invalid_uint)
for (unsigned int i=0; i != n_nodes; ++i)
elem->point(i)(2) =
libmesh_real((*elem_subsolutions[_mesh_z_var])(i));
// }
// FIXME - If the coordinate data is not in our own system, someone
// had better get around to implementing that... - RHS
// else
// {
// libmesh_not_implemented();
// }
}
Reimplemented from DiffContext.
Definition at line 490 of file fem_context.C.
References _mesh_sys, elem_fe_reinit(), elem_position_set(), and libMesh::invalid_uint.
{
// Handle a moving element if necessary
if (_mesh_sys != libMesh::invalid_uint)
{
// FIXME - not threadsafe yet!
elem_position_set(theta);
elem_fe_reinit();
}
}
Definition at line 527 of file fem_context.C.
References elem, side, and side_fe.
Referenced by FEMSystem::assembly(), and elem_side_reinit().
{
// Initialize all the interior FE objects on elem/side.
// Logging of FE::reinit is done in the FE functions
std::map<FEType, FEBase *>::iterator fe_end = side_fe.end();
for (std::map<FEType, FEBase *>::iterator i = side_fe.begin();
i != fe_end; ++i)
{
i->second->reinit(elem, side);
}
}
Reimplemented from DiffContext.
Definition at line 502 of file fem_context.C.
References _mesh_sys, elem_position_set(), elem_side_fe_reinit(), and libMesh::invalid_uint.
{
// Handle a moving element if necessary
if (_mesh_sys != libMesh::invalid_uint)
{
// FIXME - not threadsafe yet!
elem_position_set(theta);
elem_side_fe_reinit();
}
}
Definition at line 330 of file fem_context.C.
References TypeVector< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_fixed_subsolutions, and element_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealGradient> > &dphi =
element_fe_var[var]->get_dphi();
// Accumulate solution derivatives
Gradient du;
for (unsigned int l=0; l != n_dofs; l++)
du.add_scaled(dphi[l][qp], coef(l));
return du;
}
Definition at line 357 of file fem_context.C.
References TypeTensor< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_fixed_subsolutions, and element_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealTensor> > &d2phi =
element_fe_var[var]->get_d2phi();
// Accumulate solution second derivatives
Tensor d2u;
for (unsigned int l=0; l != n_dofs; l++)
d2u.add_scaled(d2phi[l][qp], coef(l));
return d2u;
}
Definition at line 304 of file fem_context.C.
References DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_fixed_subsolutions, and element_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<Real> > &phi =
element_fe_var[var]->get_phi();
// Accumulate solution value
Number u = 0.;
for (unsigned int l=0; l != n_dofs; l++)
u += phi[l][qp] * coef(l);
return u;
}
Definition at line 464 of file fem_context.C.
References dim, DiffContext::dof_indices, DiffContext::dof_indices_var, elem, DiffContext::elem_fixed_subsolutions, element_fe_var, FEInterface::inverse_map(), and FEInterface::shape().
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
Number u = 0.;
FEType fe_type = element_fe_var[var]->get_fe_type();
Point p_master = FEInterface::inverse_map(dim, fe_type, elem, p);
for (unsigned int l=0; l != n_dofs; l++)
u += FEInterface::shape(dim, fe_type, elem, l, p_master)
* coef(l);
return u;
}
Definition at line 410 of file fem_context.C.
References TypeVector< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_fixed_subsolutions, and side_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealGradient> > &dphi =
side_fe_var[var]->get_dphi();
// Accumulate solution derivatives
Gradient du;
for (unsigned int l=0; l != n_dofs; l++)
du.add_scaled(dphi[l][qp], coef(l));
return du;
}
Definition at line 437 of file fem_context.C.
References TypeTensor< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_fixed_subsolutions, and side_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealTensor> > &d2phi =
side_fe_var[var]->get_d2phi();
// Accumulate solution second derivatives
Tensor d2u;
for (unsigned int l=0; l != n_dofs; l++)
d2u.add_scaled(d2phi[l][qp], coef(l));
return d2u;
}
Definition at line 384 of file fem_context.C.
References DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_fixed_subsolutions, and side_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_fixed_subsolutions.size() > var);
libmesh_assert (elem_fixed_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_fixed_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<Real> > &phi =
side_fe_var[var]->get_phi();
// Accumulate solution value
Number u = 0.;
for (unsigned int l=0; l != n_dofs; l++)
u += phi[l][qp] * coef(l);
return u;
}
Definition at line 144 of file fem_context.C.
References TypeVector< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_subsolutions, and element_fe_var.
Referenced by FEMSystem::eulerian_residual().
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealGradient> > &dphi =
element_fe_var[var]->get_dphi();
// Accumulate solution derivatives
Gradient du;
for (unsigned int l=0; l != n_dofs; l++)
du.add_scaled(dphi[l][qp], coef(l));
return du;
}
Definition at line 171 of file fem_context.C.
References TypeTensor< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_subsolutions, and element_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealTensor> > &d2phi =
element_fe_var[var]->get_d2phi();
// Accumulate solution second derivatives
Tensor d2u;
for (unsigned int l=0; l != n_dofs; l++)
d2u.add_scaled(d2phi[l][qp], coef(l));
return d2u;
}
Definition at line 118 of file fem_context.C.
References DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_subsolutions, and element_fe_var.
Referenced by FEMSystem::mass_residual().
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<Real> > &phi =
element_fe_var[var]->get_phi();
// Accumulate solution value
Number u = 0.;
for (unsigned int l=0; l != n_dofs; l++)
u += phi[l][qp] * coef(l);
return u;
}
Definition at line 278 of file fem_context.C.
References dim, DiffContext::dof_indices, DiffContext::dof_indices_var, elem, DiffContext::elem_subsolutions, element_fe_var, FEInterface::inverse_map(), and FEInterface::shape().
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
Number u = 0.;
FEType fe_type = element_fe_var[var]->get_fe_type();
Point p_master = FEInterface::inverse_map(dim, fe_type, elem, p);
for (unsigned int l=0; l != n_dofs; l++)
u += FEInterface::shape(dim, fe_type, elem, l, p_master)
* coef(l);
return u;
}
Definition at line 653 of file fem_context.C.
References _mesh_sys, System::current_solution(), DiffContext::dof_indices, DofMap::dof_indices(), DiffContext::dof_indices_var, elem, elem_fe_reinit(), DiffContext::elem_fixed_solution, DiffContext::elem_fixed_subsolutions, DiffContext::elem_jacobian, elem_position_set(), DiffContext::elem_residual, DiffContext::elem_solution, DiffContext::elem_subjacobians, DiffContext::elem_subresiduals, DiffContext::elem_subsolutions, System::get_dof_map(), libMesh::invalid_uint, System::n_vars(), DenseMatrix< T >::resize(), DenseVector< T >::resize(), and DifferentiableSystem::use_fixed_solution.
Referenced by FEMSystem::assemble_qoi(), FEMSystem::assemble_qoi_derivative(), FEMSystem::assembly(), FEMSystem::mesh_position_set(), and FEMSystem::postprocess().
{
elem = e;
// Initialize the per-element data for elem.
sys.get_dof_map().dof_indices (elem, dof_indices);
unsigned int n_dofs = dof_indices.size();
elem_solution.resize(n_dofs);
if (sys.use_fixed_solution)
elem_fixed_solution.resize(n_dofs);
for (unsigned int i=0; i != n_dofs; ++i)
elem_solution(i) = sys.current_solution(dof_indices[i]);
// These resize calls also zero out the residual and jacobian
elem_residual.resize(n_dofs);
elem_jacobian.resize(n_dofs, n_dofs);
// Initialize the per-variable data for elem.
unsigned int sub_dofs = 0;
for (unsigned int i=0; i != sys.n_vars(); ++i)
{
sys.get_dof_map().dof_indices (elem, dof_indices_var[i], i);
elem_subsolutions[i]->reposition
(sub_dofs, dof_indices_var[i].size());
if (sys.use_fixed_solution)
elem_fixed_subsolutions[i]->reposition
(sub_dofs, dof_indices_var[i].size());
elem_subresiduals[i]->reposition
(sub_dofs, dof_indices_var[i].size());
for (unsigned int j=0; j != i; ++j)
{
elem_subjacobians[i][j]->reposition
(sub_dofs, elem_subresiduals[j]->i_off(),
dof_indices_var[i].size(),
dof_indices_var[j].size());
elem_subjacobians[j][i]->reposition
(elem_subresiduals[j]->i_off(), sub_dofs,
dof_indices_var[j].size(),
dof_indices_var[i].size());
}
elem_subjacobians[i][i]->reposition
(sub_dofs, sub_dofs,
dof_indices_var[i].size(),
dof_indices_var[i].size());
sub_dofs += dof_indices_var[i].size();
}
libmesh_assert(sub_dofs == n_dofs);
// Moving the mesh may be necessary
if (_mesh_sys != libMesh::invalid_uint)
elem_position_set(1.);
// Reinitializing the FE objects is definitely necessary
this->elem_fe_reinit();
}
Definition at line 224 of file fem_context.C.
References TypeVector< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_subsolutions, and side_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealGradient> > &dphi =
side_fe_var[var]->get_dphi();
// Accumulate solution derivatives
Gradient du;
for (unsigned int l=0; l != n_dofs; l++)
du.add_scaled(dphi[l][qp], coef(l));
return du;
}
Definition at line 251 of file fem_context.C.
References TypeTensor< T >::add_scaled(), DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_subsolutions, and side_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<RealTensor> > &d2phi =
side_fe_var[var]->get_d2phi();
// Accumulate solution second derivatives
Tensor d2u;
for (unsigned int l=0; l != n_dofs; l++)
d2u.add_scaled(d2phi[l][qp], coef(l));
return d2u;
}
Definition at line 198 of file fem_context.C.
References DiffContext::dof_indices, DiffContext::dof_indices_var, DiffContext::elem_subsolutions, and side_fe_var.
{
// Get local-to-global dof index lookup
libmesh_assert (dof_indices.size() > var);
const unsigned int n_dofs = dof_indices_var[var].size();
// Get current local coefficients
libmesh_assert (elem_subsolutions.size() > var);
libmesh_assert (elem_subsolutions[var] != NULL);
DenseSubVector<Number> &coef = *elem_subsolutions[var];
// Get shape function values at quadrature point
const std::vector<std::vector<Real> > &phi =
side_fe_var[var]->get_phi();
// Accumulate solution value
Number u = 0.;
for (unsigned int l=0; l != n_dofs; l++)
u += phi[l][qp] * coef(l);
return u;
}
Definition at line 231 of file fem_context.h.
Referenced by elem_position_get(), elem_position_set(), elem_reinit(), elem_side_reinit(), and reinit().
Definition at line 231 of file fem_context.h.
Referenced by elem_position_get(), and elem_position_set().
Definition at line 231 of file fem_context.h.
Referenced by elem_position_get(), and elem_position_set().
Definition at line 231 of file fem_context.h.
Referenced by elem_position_get(), and elem_position_set().
Definition at line 246 of file fem_context.h.
Referenced by FEMContext(), fixed_point_value(), and point_value().
Definition at line 148 of file diff_context.h.
Referenced by FEMSystem::assemble_qoi_derivative(), FEMSystem::assembly(), EulerSolver::element_residual(), Euler2Solver::element_residual(), fixed_interior_gradient(), fixed_interior_hessian(), fixed_interior_value(), fixed_point_value(), fixed_side_gradient(), fixed_side_hessian(), fixed_side_value(), interior_gradient(), interior_hessian(), interior_value(), FEMSystem::numerical_jacobian(), point_value(), reinit(), side_gradient(), side_hessian(), EulerSolver::side_residual(), Euler2Solver::side_residual(), and side_value().
Definition at line 149 of file diff_context.h.
Referenced by FEMSystem::eulerian_residual(), fixed_interior_gradient(), fixed_interior_hessian(), fixed_interior_value(), fixed_point_value(), fixed_side_gradient(), fixed_side_hessian(), fixed_side_value(), interior_gradient(), interior_hessian(), interior_value(), FEMSystem::mass_residual(), FEMSystem::mesh_position_get(), FEMSystem::numerical_jacobian(), point_value(), reinit(), side_gradient(), side_hessian(), and side_value().
Definition at line 236 of file fem_context.h.
Referenced by FEMSystem::assemble_qoi(), FEMSystem::assemble_qoi_derivative(), FEMSystem::assembly(), elem_fe_reinit(), elem_position_get(), elem_position_set(), elem_side_fe_reinit(), fixed_point_value(), FEMSystem::mesh_position_get(), FEMSystem::mesh_position_set(), FEMSystem::numerical_jacobian(), point_value(), FEMSystem::postprocess(), and reinit().
Definition at line 105 of file diff_context.h.
Referenced by DiffContext::DiffContext(), SteadySolver::element_residual(), EulerSolver::element_residual(), Euler2Solver::element_residual(), reinit(), SteadySolver::side_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
Definition at line 106 of file diff_context.h.
Referenced by DiffContext::DiffContext(), fixed_interior_gradient(), fixed_interior_hessian(), fixed_interior_value(), fixed_point_value(), fixed_side_gradient(), fixed_side_hessian(), fixed_side_value(), reinit(), and DiffContext::~DiffContext().
Definition at line 137 of file diff_context.h.
Referenced by FEMSystem::assembly(), DiffContext::DiffContext(), EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::numerical_jacobian(), reinit(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
Definition at line 131 of file diff_context.h.
Referenced by FEMSystem::assemble_qoi_derivative(), FEMSystem::assembly(), DiffContext::DiffContext(), EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::numerical_jacobian(), reinit(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
Definition at line 97 of file diff_context.h.
Referenced by DiffContext::DiffContext(), SteadySolver::element_residual(), EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::numerical_jacobian(), reinit(), SteadySolver::side_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
Definition at line 119 of file diff_context.h.
Referenced by EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::mass_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
Definition at line 143 of file diff_context.h.
Referenced by DiffContext::DiffContext(), FEMSystem::eulerian_residual(), FEMSystem::mass_residual(), reinit(), and DiffContext::~DiffContext().
Definition at line 142 of file diff_context.h.
Referenced by DiffContext::DiffContext(), FEMSystem::eulerian_residual(), FEMSystem::mass_residual(), reinit(), and DiffContext::~DiffContext().
Definition at line 98 of file diff_context.h.
Referenced by DiffContext::DiffContext(), elem_position_get(), elem_position_set(), interior_gradient(), interior_hessian(), interior_value(), FEMSystem::mesh_position_get(), point_value(), reinit(), side_gradient(), side_hessian(), side_value(), and DiffContext::~DiffContext().
Definition at line 197 of file fem_context.h.
Referenced by elem_fe_reinit(), FEMContext(), and ~FEMContext().
Definition at line 204 of file fem_context.h.
Referenced by elem_position_get(), elem_position_set(), FEMSystem::eulerian_residual(), FEMContext(), fixed_interior_gradient(), fixed_interior_hessian(), fixed_interior_value(), fixed_point_value(), FEMSystem::init_context(), interior_gradient(), interior_hessian(), interior_value(), FEMSystem::mass_residual(), and point_value().
Definition at line 212 of file fem_context.h.
Referenced by FEMSystem::eulerian_residual(), FEMContext(), FEMSystem::mass_residual(), and ~FEMContext().
Definition at line 126 of file diff_context.h.
Referenced by SteadySolver::element_residual(), EulerSolver::element_residual(), Euler2Solver::element_residual(), SteadySolver::side_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
Definition at line 71 of file diff_context.h.
Definition at line 112 of file diff_context.h.
Definition at line 241 of file fem_context.h.
Referenced by FEMSystem::assemble_qoi(), FEMSystem::assemble_qoi_derivative(), FEMSystem::assembly(), elem_side_fe_reinit(), and FEMSystem::postprocess().
Definition at line 198 of file fem_context.h.
Referenced by FEMSystem::assemble_qoi(), FEMSystem::assemble_qoi_derivative(), elem_side_fe_reinit(), FEMContext(), FEMSystem::postprocess(), and ~FEMContext().
Definition at line 205 of file fem_context.h.
Referenced by FEMContext(), fixed_side_gradient(), fixed_side_hessian(), fixed_side_value(), side_gradient(), side_hessian(), and side_value().
Definition at line 213 of file fem_context.h.
Referenced by FEMContext(), and ~FEMContext().
Definition at line 91 of file diff_context.h.
Generated automatically by Doxygen for libMesh from the source code.