]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement FE::convert_generalized_support_point_values_to_nodal_values() for FE_Q_iso_Q1. 4202/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 6 Apr 2017 14:52:31 +0000 (08:52 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 6 Apr 2017 14:52:31 +0000 (08:52 -0600)
include/deal.II/fe/fe_q_iso_q1.h
source/fe/fe_q_iso_q1.cc

index d163801c97b639f80ff39de816a9f0f9b67a1c53..92bf7f5cb82fa22a2c5f20a2a0ff2eace5ac853a 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.
 //
@@ -121,6 +121,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;
+
   /**
    * @name Functions to support hp
    * @{
index 6e362f651869bc7b1c1fe3bb54b91fa904d19541..fdbf224aa57128c6260e5b5cef6da0a7d8fd6636 100644 (file)
@@ -15,6 +15,7 @@
 
 
 #include <deal.II/base/quadrature_lib.h>
+#include <deal.II/lac/vector.h>
 #include <deal.II/fe/fe_q_iso_q1.h>
 #include <deal.II/fe/fe_nothing.h>
 
@@ -67,6 +68,29 @@ FE_Q_iso_Q1<dim,spacedim>::get_name () const
 
 
 
+template <int dim, int spacedim>
+void
+FE_Q_iso_Q1<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_iso_Q1<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.