]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement FE::convert_generalized_support_point_values_to_nodal_values() for FE_Q... 4185/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 5 Apr 2017 16:38:11 +0000 (10:38 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Wed, 5 Apr 2017 21:32:09 +0000 (15:32 -0600)
doc/news/changes/minor/20170405Bangerth [new file with mode: 0644]
include/deal.II/fe/fe_dgq.h
include/deal.II/fe/fe_q.h
include/deal.II/fe/fe_q_base.h
source/fe/fe_dgq.cc
source/fe/fe_q.cc

diff --git a/doc/news/changes/minor/20170405Bangerth b/doc/news/changes/minor/20170405Bangerth
new file mode 100644 (file)
index 0000000..d9c9c34
--- /dev/null
@@ -0,0 +1,6 @@
+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)
index 1c04563d412d9f3095b20d71cad75bd6f9d05c14..693cc52bca8c8894280b15490cc1b3fe2957277d 100644 (file)
@@ -300,6 +300,18 @@ public:
   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.
@@ -409,6 +421,18 @@ public:
    */
   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.
index a19aa9bf4ea2f9f224e8cafe8b94ab6574968872..7c08e61e45b71e48059fa106fb16b39351e55542 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// 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.
 //
@@ -586,6 +586,18 @@ public:
    */
   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:
 
   /**
index 2a39e81e90c5775f501a3a594dd88f69bbcb76eb..4b9e35e82e6ff2e5ca355a04fe221c65c72e153a 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// 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.
 //
index af344ba3fb1dc315637b54a37aab4026e0e2fa44..bc9159982d4ec14d05eb400e0d5c30316ae7f896 100644 (file)
@@ -16,6 +16,7 @@
 
 #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>
@@ -102,6 +103,29 @@ FE_DGQ<dim, spacedim>::get_name () const
 
 
 
+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
@@ -759,6 +783,30 @@ FE_DGQArbitraryNodes<dim,spacedim>::get_name () 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
index af64fe48c6f8505fcd0258e433d39279999a8615..f0e272f4a4f6fe46d8f0fc97f98322a6434decce 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// 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.
 //
@@ -15,6 +15,7 @@
 
 
 #include <deal.II/base/quadrature_lib.h>
+#include <deal.II/lac/vector.h>
 #include <deal.II/fe/fe_q.h>
 
 #include <vector>
@@ -132,6 +133,29 @@ FE_Q<dim,spacedim>::get_name () const
 
 
 
+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

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.