--- /dev/null
+New: The classes FE_Q, FE_DGQ, and FE_DGQArbitraryNodes now implement
+the
+FiniteElement::convert_generalized_support_point_values_to_nodal_values()
+interface.
+<br>
+(Wolfgang Bangerth, 2017/04/05)
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
+ /**
+ * Implementation of the corresponding function in the FiniteElement
+ * class. Since the current element is interpolatory, the nodal
+ * values are exactly the support point values. Furthermore, since
+ * the current element is scalar, the support point values need to
+ * be vectors of length 1.
+ */
+ virtual
+ void
+ convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
+ std::vector<double> &nodal_values) const;
+
/**
* Determine an estimate for the memory consumption (in bytes) of this
* object.
*/
virtual std::string get_name () const;
+ /**
+ * Implementation of the corresponding function in the FiniteElement
+ * class. Since the current element is interpolatory, the nodal
+ * values are exactly the support point values. Furthermore, since
+ * the current element is scalar, the support point values need to
+ * be vectors of length 1.
+ */
+ virtual
+ void
+ convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
+ std::vector<double> &nodal_values) const;
+
protected:
/**
* @p clone function instead of a copy constructor.
// ---------------------------------------------------------------------
//
-// Copyright (C) 2000 - 2016 by the deal.II authors
+// Copyright (C) 2000 - 2017 by the deal.II authors
//
// This file is part of the deal.II library.
//
*/
virtual std::string get_name () const;
+ /**
+ * Implementation of the corresponding function in the FiniteElement
+ * class. Since the current element is interpolatory, the nodal
+ * values are exactly the support point values. Furthermore, since
+ * the current element is scalar, the support point values need to
+ * be vectors of length 1.
+ */
+ virtual
+ void
+ convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
+ std::vector<double> &nodal_values) const;
+
protected:
/**
// ---------------------------------------------------------------------
//
-// Copyright (C) 2000 - 2016 by the deal.II authors
+// Copyright (C) 2000 - 2017 by the deal.II authors
//
// This file is part of the deal.II library.
//
#include <deal.II/base/quadrature.h>
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/lac/vector.h>
#include <deal.II/fe/fe.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_tools.h>
+template <int dim, int spacedim>
+void
+FE_DGQ<dim,spacedim>::
+convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
+ std::vector<double> &nodal_values) const
+{
+ AssertDimension (support_point_values.size(),
+ this->get_unit_support_points().size());
+ AssertDimension (support_point_values.size(),
+ nodal_values.size());
+ AssertDimension (this->dofs_per_cell,
+ nodal_values.size());
+
+ for (unsigned int i=0; i<this->dofs_per_cell; ++i)
+ {
+ AssertDimension (support_point_values[i].size(), 1);
+
+ nodal_values[i] = support_point_values[i](0);
+ }
+}
+
+
+
template <int dim, int spacedim>
FiniteElement<dim,spacedim> *
FE_DGQ<dim, spacedim>::clone() const
+template <int dim, int spacedim>
+void
+FE_DGQArbitraryNodes<dim,spacedim>::
+convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
+ std::vector<double> &nodal_values) const
+{
+ AssertDimension (support_point_values.size(),
+ this->get_unit_support_points().size());
+ AssertDimension (support_point_values.size(),
+ nodal_values.size());
+ AssertDimension (this->dofs_per_cell,
+ nodal_values.size());
+
+ for (unsigned int i=0; i<this->dofs_per_cell; ++i)
+ {
+ AssertDimension (support_point_values[i].size(), 1);
+
+ nodal_values[i] = support_point_values[i](0);
+ }
+}
+
+
+
+
template <int dim, int spacedim>
FiniteElement<dim,spacedim> *
FE_DGQArbitraryNodes<dim,spacedim>::clone() const
// ---------------------------------------------------------------------
//
-// Copyright (C) 2000 - 2016 by the deal.II authors
+// Copyright (C) 2000 - 2017 by the deal.II authors
//
// This file is part of the deal.II library.
//
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/lac/vector.h>
#include <deal.II/fe/fe_q.h>
#include <vector>
+template <int dim, int spacedim>
+void
+FE_Q<dim,spacedim>::
+convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
+ std::vector<double> &nodal_values) const
+{
+ AssertDimension (support_point_values.size(),
+ this->get_unit_support_points().size());
+ AssertDimension (support_point_values.size(),
+ nodal_values.size());
+ AssertDimension (this->dofs_per_cell,
+ nodal_values.size());
+
+ for (unsigned int i=0; i<this->dofs_per_cell; ++i)
+ {
+ AssertDimension (support_point_values[i].size(), 1);
+
+ nodal_values[i] = support_point_values[i](0);
+ }
+}
+
+
+
template <int dim, int spacedim>
FiniteElement<dim,spacedim> *
FE_Q<dim,spacedim>::clone() const