*
* @author Joshua White, Wolfgang Bangerth
*/
-template <int dim>
-class FE_Nothing : public FiniteElement<dim>
+template <int dim, int spacedim=dim>
+class FE_Nothing : public FiniteElement<dim,spacedim>
{
public:
* number of components to give this
* finite element (default = 1).
*/
- FE_Nothing (unsigned int n_components = 1);
+ FE_Nothing (const unsigned int n_components = 1);
/**
* A sort of virtual copy
* do so through this function.
*/
virtual
- FiniteElement<dim> *
+ FiniteElement<dim,spacedim> *
clone() const;
/**
*/
virtual
void
- fill_fe_values (const Mapping<dim> &mapping,
- const typename Triangulation<dim>::cell_iterator &cell,
+ fill_fe_values (const Mapping<dim,spacedim> &mapping,
+ const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
- typename Mapping<dim>::InternalDataBase &mapping_data,
- typename Mapping<dim>::InternalDataBase &fedata,
- FEValuesData<dim,dim> &data,
+ typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+ typename Mapping<dim,spacedim>::InternalDataBase &fedata,
+ FEValuesData<dim,spacedim> &data,
CellSimilarity::Similarity &cell_similarity) const;
/**
*/
virtual
void
- fill_fe_face_values (const Mapping<dim> &mapping,
- const typename Triangulation<dim> :: cell_iterator &cell,
+ fill_fe_face_values (const Mapping<dim,spacedim> &mapping,
+ const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face,
const Quadrature<dim-1> & quadrature,
- typename Mapping<dim> :: InternalDataBase &mapping_data,
- typename Mapping<dim> :: InternalDataBase &fedata,
- FEValuesData<dim,dim> &data) const;
+ typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+ typename Mapping<dim,spacedim>::InternalDataBase &fedata,
+ FEValuesData<dim,spacedim> &data) const;
/**
* Fill the fields of
*/
virtual
void
- fill_fe_subface_values (const Mapping<dim> &mapping,
- const typename Triangulation<dim>::cell_iterator &cell,
+ fill_fe_subface_values (const Mapping<dim,spacedim> &mapping,
+ const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face,
const unsigned int subface,
const Quadrature<dim-1> & quadrature,
- typename Mapping<dim>::InternalDataBase &mapping_data,
- typename Mapping<dim>::InternalDataBase &fedata,
- FEValuesData<dim,dim> &data) const;
+ typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+ typename Mapping<dim,spacedim>::InternalDataBase &fedata,
+ FEValuesData<dim,spacedim> &data) const;
/**
* Prepare internal data
* element.
*/
virtual
- typename Mapping<dim>::InternalDataBase *
+ typename Mapping<dim,spacedim>::InternalDataBase *
get_data (const UpdateFlags update_flags,
- const Mapping<dim> &mapping,
+ const Mapping<dim,spacedim> &mapping,
const Quadrature<dim> &quadrature) const;
/**
*/
virtual
FiniteElementDomination::Domination
- compare_for_face_domination (const FiniteElement<dim> &fe_other) const;
+ compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const;
virtual
std::vector<std::pair<unsigned int, unsigned int> >
- hp_vertex_dof_identities (const FiniteElement<dim> &fe_other) const;
+ hp_vertex_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
virtual
std::vector<std::pair<unsigned int, unsigned int> >
- hp_line_dof_identities (const FiniteElement<dim> &fe_other) const;
+ hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
virtual
std::vector<std::pair<unsigned int, unsigned int> >
- hp_quad_dof_identities (const FiniteElement<dim> &fe_other) const;
+ hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
virtual
bool
virtual
void
- get_face_interpolation_matrix (const FiniteElement<dim> &source_fe,
+ get_face_interpolation_matrix (const FiniteElement<dim,spacedim> &source_fe,
FullMatrix<double> &interpolation_matrix) const;
virtual
void
- get_subface_interpolation_matrix (const FiniteElement<dim> &source_fe,
+ get_subface_interpolation_matrix (const FiniteElement<dim,spacedim> &source_fe,
const unsigned int index,
FullMatrix<double> &interpolation_matrix) const;
-template <int dim>
-FE_Nothing<dim>::FE_Nothing (const unsigned int n_components)
+template <int dim, int spacedim>
+FE_Nothing<dim,spacedim>::FE_Nothing (const unsigned int n_components)
:
- FiniteElement<dim>
+ FiniteElement<dim,spacedim>
(FiniteElementData<dim>(std::vector<unsigned>(dim+1,0),
n_components, 0,
FiniteElementData<dim>::unknown),
}
-template <int dim>
-FiniteElement<dim> *
-FE_Nothing<dim>::clone() const
+template <int dim, int spacedim>
+FiniteElement<dim,spacedim> *
+FE_Nothing<dim,spacedim>::clone() const
{
- return new FE_Nothing<dim>(*this);
+ return new FE_Nothing<dim,spacedim>(*this);
}
-template <int dim>
+template <int dim, int spacedim>
std::string
-FE_Nothing<dim>::get_name () const
+FE_Nothing<dim,spacedim>::get_name () const
{
std::ostringstream namebuf;
namebuf << "FE_Nothing<" << dim << ">(";
-template <int dim>
+template <int dim, int spacedim>
UpdateFlags
-FE_Nothing<dim>::update_once (const UpdateFlags /*flags*/) const
+FE_Nothing<dim,spacedim>::update_once (const UpdateFlags /*flags*/) const
{
return update_default;
}
-template <int dim>
+template <int dim, int spacedim>
UpdateFlags
-FE_Nothing<dim>::update_each (const UpdateFlags /*flags*/) const
+FE_Nothing<dim,spacedim>::update_each (const UpdateFlags /*flags*/) const
{
return update_default;
}
-template <int dim>
+template <int dim, int spacedim>
double
-FE_Nothing<dim>::shape_value (const unsigned int /*i*/,
- const Point<dim> & /*p*/) const
+FE_Nothing<dim,spacedim>::shape_value (const unsigned int /*i*/,
+ const Point<dim> & /*p*/) const
{
Assert(false,ExcMessage(zero_dof_message));
return 0;
-template <int dim>
-typename Mapping<dim>::InternalDataBase *
-FE_Nothing<dim>::get_data (const UpdateFlags /*flags*/,
- const Mapping<dim> & /*mapping*/,
- const Quadrature<dim> & /*quadrature*/) const
+template <int dim, int spacedim>
+typename Mapping<dim,spacedim>::InternalDataBase *
+FE_Nothing<dim,spacedim>::get_data (const UpdateFlags /*flags*/,
+ const Mapping<dim,spacedim> & /*mapping*/,
+ const Quadrature<dim> & /*quadrature*/) const
{
// Create a default data object. Normally we would then
// need to resize things to hold the appropriate numbers
// of dofs, but in this case all data fields are empty.
- typename Mapping<dim>::InternalDataBase *data
- = new typename FiniteElement<dim>::InternalDataBase();
+ typename Mapping<dim,spacedim>::InternalDataBase *data
+ = new typename FiniteElement<dim,spacedim>::InternalDataBase();
return data;
}
-template <int dim>
+template <int dim, int spacedim>
void
-FE_Nothing<dim>::
-fill_fe_values (const Mapping<dim> & /*mapping*/,
- const typename Triangulation<dim>::cell_iterator & /*cell*/,
+FE_Nothing<dim,spacedim>::
+fill_fe_values (const Mapping<dim,spacedim> & /*mapping*/,
+ const typename Triangulation<dim,spacedim>::cell_iterator & /*cell*/,
const Quadrature<dim> & /*quadrature*/,
- typename Mapping<dim>::InternalDataBase & /*mapping_data*/,
- typename Mapping<dim>::InternalDataBase & /*fedata*/,
- FEValuesData<dim,dim> & /*data*/,
+ typename Mapping<dim,spacedim>::InternalDataBase & /*mapping_data*/,
+ typename Mapping<dim,spacedim>::InternalDataBase & /*fedata*/,
+ FEValuesData<dim,spacedim> & /*data*/,
CellSimilarity::Similarity & /*cell_similarity*/) const
{
// leave data fields empty
-template <int dim>
+template <int dim, int spacedim>
void
-FE_Nothing<dim>::
-fill_fe_face_values (const Mapping<dim> & /*mapping*/,
- const typename Triangulation<dim>::cell_iterator & /*cell*/,
+FE_Nothing<dim,spacedim>::
+fill_fe_face_values (const Mapping<dim,spacedim> & /*mapping*/,
+ const typename Triangulation<dim,spacedim>::cell_iterator & /*cell*/,
const unsigned int /*face*/,
const Quadrature<dim-1> & /*quadrature*/,
- typename Mapping<dim>::InternalDataBase & /*mapping_data*/,
- typename Mapping<dim>::InternalDataBase & /*fedata*/,
- FEValuesData<dim,dim> & /*data*/) const
+ typename Mapping<dim,spacedim>::InternalDataBase & /*mapping_data*/,
+ typename Mapping<dim,spacedim>::InternalDataBase & /*fedata*/,
+ FEValuesData<dim,spacedim> & /*data*/) const
{
// leave data fields empty
}
-template <int dim>
+template <int dim, int spacedim>
void
-FE_Nothing<dim>::
-fill_fe_subface_values (const Mapping<dim> & /*mapping*/,
- const typename Triangulation<dim>::cell_iterator & /*cell*/,
+FE_Nothing<dim,spacedim>::
+fill_fe_subface_values (const Mapping<dim,spacedim> & /*mapping*/,
+ const typename Triangulation<dim,spacedim>::cell_iterator & /*cell*/,
const unsigned int /*face*/,
const unsigned int /*subface*/,
const Quadrature<dim-1> & /*quadrature*/,
- typename Mapping<dim>::InternalDataBase & /*mapping_data*/,
- typename Mapping<dim>::InternalDataBase & /*fedata*/,
- FEValuesData<dim,dim> & /*data*/) const
+ typename Mapping<dim,spacedim>::InternalDataBase & /*mapping_data*/,
+ typename Mapping<dim,spacedim>::InternalDataBase & /*fedata*/,
+ FEValuesData<dim,spacedim> & /*data*/) const
{
// leave data fields empty
}
-template <int dim>
+template <int dim, int spacedim>
FiniteElementDomination::Domination
-FE_Nothing<dim> ::
-compare_for_face_domination (const FiniteElement<dim> &) const
+FE_Nothing<dim,spacedim> ::
+compare_for_face_domination (const FiniteElement<dim,spacedim> &) const
{
return FiniteElementDomination::no_requirements;
}
-template <int dim>
+template <int dim, int spacedim>
std::vector<std::pair<unsigned int, unsigned int> >
-FE_Nothing<dim> ::
-hp_vertex_dof_identities (const FiniteElement<dim> &/*fe_other*/) const
+FE_Nothing<dim,spacedim> ::
+hp_vertex_dof_identities (const FiniteElement<dim,spacedim> &/*fe_other*/) const
{
// the FE_Nothing has no
// degrees of freedom, so there
}
-template <int dim>
+template <int dim, int spacedim>
std::vector<std::pair<unsigned int, unsigned int> >
-FE_Nothing<dim> ::
-hp_line_dof_identities (const FiniteElement<dim> &/*fe_other*/) const
+FE_Nothing<dim,spacedim> ::
+hp_line_dof_identities (const FiniteElement<dim,spacedim> &/*fe_other*/) const
{
// the FE_Nothing has no
// degrees of freedom, so there
}
-template <int dim>
+template <int dim, int spacedim>
std::vector<std::pair<unsigned int, unsigned int> >
-FE_Nothing<dim> ::
-hp_quad_dof_identities (const FiniteElement<dim> &/*fe_other*/) const
+FE_Nothing<dim,spacedim> ::
+hp_quad_dof_identities (const FiniteElement<dim,spacedim> &/*fe_other*/) const
{
// the FE_Nothing has no
// degrees of freedom, so there
}
-template <int dim>
+template <int dim, int spacedim>
bool
-FE_Nothing<dim> ::
+FE_Nothing<dim,spacedim> ::
hp_constraints_are_implemented () const
{
return true;
}
-template <int dim>
+template <int dim, int spacedim>
void
-FE_Nothing<dim>::
-get_face_interpolation_matrix (const FiniteElement<dim> &/*source_fe*/,
+FE_Nothing<dim,spacedim>::
+get_face_interpolation_matrix (const FiniteElement<dim,spacedim> &/*source_fe*/,
FullMatrix<double> &interpolation_matrix) const
{
// since this element has no face dofs, the
}
-template <int dim>
+template <int dim, int spacedim>
void
-FE_Nothing<dim>::
-get_subface_interpolation_matrix (const FiniteElement<dim> & /*source_fe*/,
+FE_Nothing<dim,spacedim>::
+get_subface_interpolation_matrix (const FiniteElement<dim,spacedim> & /*source_fe*/,
const unsigned int /*index*/,
FullMatrix<double> &interpolation_matrix) const
{