From e31595757f587a2bceb4928d4ac8147852c0a60d Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Wed, 5 Nov 2014 11:09:43 -0600 Subject: [PATCH] Add second template parameter to FE_Nothing. For historical (and probably not very good) reasons, finite element classes have a second template argument denoting the space dimension. FE_Nothing never got this argument and so it wasn't possible to use it in codimension-one situations. This patch adds the second argument to make it possible to use it in these situations. There is not much that one can test this element with, other than ensure that it compiles. The next commit will actually use the element in a concrete situation. --- doc/news/changes.h | 8 ++ include/deal.II/fe/fe_nothing.h | 54 +++++++------- source/fe/fe_nothing.cc | 128 ++++++++++++++++---------------- source/fe/fe_nothing.inst.in | 6 +- 4 files changed, 103 insertions(+), 93 deletions(-) diff --git a/doc/news/changes.h b/doc/news/changes.h index cc7b71d160..8bc18b2d41 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -333,6 +333,14 @@ inconvenience this causes.

Specific improvements

    +
  1. New: The FE_Nothing class now has a second template argument + corresponding to the space dimension in which the mesh is embedded, + like many other classes. This allows to use this element in codimension + one cases as well now. +
    + (Wolfgang Bangerth, 2014/11/05) +
  2. +
  3. Fixed: Using the FEEvaluation framework did not work for scalar elements in 1d because there were conflicting partial specializations. This is now fixed. diff --git a/include/deal.II/fe/fe_nothing.h b/include/deal.II/fe/fe_nothing.h index 26a4234494..e4ef19ccd4 100644 --- a/include/deal.II/fe/fe_nothing.h +++ b/include/deal.II/fe/fe_nothing.h @@ -76,8 +76,8 @@ DEAL_II_NAMESPACE_OPEN * * @author Joshua White, Wolfgang Bangerth */ -template -class FE_Nothing : public FiniteElement +template +class FE_Nothing : public FiniteElement { public: @@ -86,7 +86,7 @@ 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 @@ -100,7 +100,7 @@ public: * do so through this function. */ virtual - FiniteElement * + FiniteElement * clone() const; /** @@ -190,12 +190,12 @@ public: */ virtual void - fill_fe_values (const Mapping &mapping, - const typename Triangulation::cell_iterator &cell, + fill_fe_values (const Mapping &mapping, + const typename Triangulation::cell_iterator &cell, const Quadrature &quadrature, - typename Mapping::InternalDataBase &mapping_data, - typename Mapping::InternalDataBase &fedata, - FEValuesData &data, + typename Mapping::InternalDataBase &mapping_data, + typename Mapping::InternalDataBase &fedata, + FEValuesData &data, CellSimilarity::Similarity &cell_similarity) const; /** @@ -212,13 +212,13 @@ public: */ virtual void - fill_fe_face_values (const Mapping &mapping, - const typename Triangulation :: cell_iterator &cell, + fill_fe_face_values (const Mapping &mapping, + const typename Triangulation::cell_iterator &cell, const unsigned int face, const Quadrature & quadrature, - typename Mapping :: InternalDataBase &mapping_data, - typename Mapping :: InternalDataBase &fedata, - FEValuesData &data) const; + typename Mapping::InternalDataBase &mapping_data, + typename Mapping::InternalDataBase &fedata, + FEValuesData &data) const; /** * Fill the fields of @@ -234,14 +234,14 @@ public: */ virtual void - fill_fe_subface_values (const Mapping &mapping, - const typename Triangulation::cell_iterator &cell, + fill_fe_subface_values (const Mapping &mapping, + const typename Triangulation::cell_iterator &cell, const unsigned int face, const unsigned int subface, const Quadrature & quadrature, - typename Mapping::InternalDataBase &mapping_data, - typename Mapping::InternalDataBase &fedata, - FEValuesData &data) const; + typename Mapping::InternalDataBase &mapping_data, + typename Mapping::InternalDataBase &fedata, + FEValuesData &data) const; /** * Prepare internal data @@ -260,9 +260,9 @@ public: * element. */ virtual - typename Mapping::InternalDataBase * + typename Mapping::InternalDataBase * get_data (const UpdateFlags update_flags, - const Mapping &mapping, + const Mapping &mapping, const Quadrature &quadrature) const; /** @@ -285,21 +285,21 @@ public: */ virtual FiniteElementDomination::Domination - compare_for_face_domination (const FiniteElement &fe_other) const; + compare_for_face_domination (const FiniteElement &fe_other) const; virtual std::vector > - hp_vertex_dof_identities (const FiniteElement &fe_other) const; + hp_vertex_dof_identities (const FiniteElement &fe_other) const; virtual std::vector > - hp_line_dof_identities (const FiniteElement &fe_other) const; + hp_line_dof_identities (const FiniteElement &fe_other) const; virtual std::vector > - hp_quad_dof_identities (const FiniteElement &fe_other) const; + hp_quad_dof_identities (const FiniteElement &fe_other) const; virtual bool @@ -321,7 +321,7 @@ public: virtual void - get_face_interpolation_matrix (const FiniteElement &source_fe, + get_face_interpolation_matrix (const FiniteElement &source_fe, FullMatrix &interpolation_matrix) const; @@ -341,7 +341,7 @@ public: virtual void - get_subface_interpolation_matrix (const FiniteElement &source_fe, + get_subface_interpolation_matrix (const FiniteElement &source_fe, const unsigned int index, FullMatrix &interpolation_matrix) const; diff --git a/source/fe/fe_nothing.cc b/source/fe/fe_nothing.cc index 0e5bcea5fb..d1e3e228c1 100644 --- a/source/fe/fe_nothing.cc +++ b/source/fe/fe_nothing.cc @@ -27,10 +27,10 @@ namespace -template -FE_Nothing::FE_Nothing (const unsigned int n_components) +template +FE_Nothing::FE_Nothing (const unsigned int n_components) : - FiniteElement + FiniteElement (FiniteElementData(std::vector(dim+1,0), n_components, 0, FiniteElementData::unknown), @@ -45,18 +45,18 @@ FE_Nothing::FE_Nothing (const unsigned int n_components) } -template -FiniteElement * -FE_Nothing::clone() const +template +FiniteElement * +FE_Nothing::clone() const { - return new FE_Nothing(*this); + return new FE_Nothing(*this); } -template +template std::string -FE_Nothing::get_name () const +FE_Nothing::get_name () const { std::ostringstream namebuf; namebuf << "FE_Nothing<" << dim << ">("; @@ -68,28 +68,28 @@ FE_Nothing::get_name () const -template +template UpdateFlags -FE_Nothing::update_once (const UpdateFlags /*flags*/) const +FE_Nothing::update_once (const UpdateFlags /*flags*/) const { return update_default; } -template +template UpdateFlags -FE_Nothing::update_each (const UpdateFlags /*flags*/) const +FE_Nothing::update_each (const UpdateFlags /*flags*/) const { return update_default; } -template +template double -FE_Nothing::shape_value (const unsigned int /*i*/, - const Point & /*p*/) const +FE_Nothing::shape_value (const unsigned int /*i*/, + const Point & /*p*/) const { Assert(false,ExcMessage(zero_dof_message)); return 0; @@ -97,31 +97,31 @@ FE_Nothing::shape_value (const unsigned int /*i*/, -template -typename Mapping::InternalDataBase * -FE_Nothing::get_data (const UpdateFlags /*flags*/, - const Mapping & /*mapping*/, - const Quadrature & /*quadrature*/) const +template +typename Mapping::InternalDataBase * +FE_Nothing::get_data (const UpdateFlags /*flags*/, + const Mapping & /*mapping*/, + const Quadrature & /*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::InternalDataBase *data - = new typename FiniteElement::InternalDataBase(); + typename Mapping::InternalDataBase *data + = new typename FiniteElement::InternalDataBase(); return data; } -template +template void -FE_Nothing:: -fill_fe_values (const Mapping & /*mapping*/, - const typename Triangulation::cell_iterator & /*cell*/, +FE_Nothing:: +fill_fe_values (const Mapping & /*mapping*/, + const typename Triangulation::cell_iterator & /*cell*/, const Quadrature & /*quadrature*/, - typename Mapping::InternalDataBase & /*mapping_data*/, - typename Mapping::InternalDataBase & /*fedata*/, - FEValuesData & /*data*/, + typename Mapping::InternalDataBase & /*mapping_data*/, + typename Mapping::InternalDataBase & /*fedata*/, + FEValuesData & /*data*/, CellSimilarity::Similarity & /*cell_similarity*/) const { // leave data fields empty @@ -129,49 +129,49 @@ fill_fe_values (const Mapping & /*mapping*/, -template +template void -FE_Nothing:: -fill_fe_face_values (const Mapping & /*mapping*/, - const typename Triangulation::cell_iterator & /*cell*/, +FE_Nothing:: +fill_fe_face_values (const Mapping & /*mapping*/, + const typename Triangulation::cell_iterator & /*cell*/, const unsigned int /*face*/, const Quadrature & /*quadrature*/, - typename Mapping::InternalDataBase & /*mapping_data*/, - typename Mapping::InternalDataBase & /*fedata*/, - FEValuesData & /*data*/) const + typename Mapping::InternalDataBase & /*mapping_data*/, + typename Mapping::InternalDataBase & /*fedata*/, + FEValuesData & /*data*/) const { // leave data fields empty } -template +template void -FE_Nothing:: -fill_fe_subface_values (const Mapping & /*mapping*/, - const typename Triangulation::cell_iterator & /*cell*/, +FE_Nothing:: +fill_fe_subface_values (const Mapping & /*mapping*/, + const typename Triangulation::cell_iterator & /*cell*/, const unsigned int /*face*/, const unsigned int /*subface*/, const Quadrature & /*quadrature*/, - typename Mapping::InternalDataBase & /*mapping_data*/, - typename Mapping::InternalDataBase & /*fedata*/, - FEValuesData & /*data*/) const + typename Mapping::InternalDataBase & /*mapping_data*/, + typename Mapping::InternalDataBase & /*fedata*/, + FEValuesData & /*data*/) const { // leave data fields empty } -template +template FiniteElementDomination::Domination -FE_Nothing :: -compare_for_face_domination (const FiniteElement &) const +FE_Nothing :: +compare_for_face_domination (const FiniteElement &) const { return FiniteElementDomination::no_requirements; } -template +template std::vector > -FE_Nothing :: -hp_vertex_dof_identities (const FiniteElement &/*fe_other*/) const +FE_Nothing :: +hp_vertex_dof_identities (const FiniteElement &/*fe_other*/) const { // the FE_Nothing has no // degrees of freedom, so there @@ -181,10 +181,10 @@ hp_vertex_dof_identities (const FiniteElement &/*fe_other*/) const } -template +template std::vector > -FE_Nothing :: -hp_line_dof_identities (const FiniteElement &/*fe_other*/) const +FE_Nothing :: +hp_line_dof_identities (const FiniteElement &/*fe_other*/) const { // the FE_Nothing has no // degrees of freedom, so there @@ -194,10 +194,10 @@ hp_line_dof_identities (const FiniteElement &/*fe_other*/) const } -template +template std::vector > -FE_Nothing :: -hp_quad_dof_identities (const FiniteElement &/*fe_other*/) const +FE_Nothing :: +hp_quad_dof_identities (const FiniteElement &/*fe_other*/) const { // the FE_Nothing has no // degrees of freedom, so there @@ -207,19 +207,19 @@ hp_quad_dof_identities (const FiniteElement &/*fe_other*/) const } -template +template bool -FE_Nothing :: +FE_Nothing :: hp_constraints_are_implemented () const { return true; } -template +template void -FE_Nothing:: -get_face_interpolation_matrix (const FiniteElement &/*source_fe*/, +FE_Nothing:: +get_face_interpolation_matrix (const FiniteElement &/*source_fe*/, FullMatrix &interpolation_matrix) const { // since this element has no face dofs, the @@ -234,10 +234,10 @@ get_face_interpolation_matrix (const FiniteElement &/*source_fe*/, } -template +template void -FE_Nothing:: -get_subface_interpolation_matrix (const FiniteElement & /*source_fe*/, +FE_Nothing:: +get_subface_interpolation_matrix (const FiniteElement & /*source_fe*/, const unsigned int /*index*/, FullMatrix &interpolation_matrix) const { diff --git a/source/fe/fe_nothing.inst.in b/source/fe/fe_nothing.inst.in index f1a670baa9..228cf903ad 100644 --- a/source/fe/fe_nothing.inst.in +++ b/source/fe/fe_nothing.inst.in @@ -15,8 +15,10 @@ -for (deal_II_dimension : DIMENSIONS) +for (deal_II_dimension, deal_II_space_dimension : DIMENSIONS) { - template class FE_Nothing; +#if deal_II_dimension <= deal_II_space_dimension + template class FE_Nothing; +#endif } -- 2.39.5