]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use FE::convert_generalized_support_point_values_to_nodal_values().
authorWolfgang Bangerth <bangerth@colostate.edu>
Tue, 14 Mar 2017 00:46:45 +0000 (18:46 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Sun, 19 Mar 2017 19:02:02 +0000 (13:02 -0600)
include/deal.II/fe/fe_tools.h
source/fe/fe_tools.cc

index e1b50ecebe5f4a153b2b06ef4b41cc8c52534f32..b4f8a0270bfa1a3cac07490e4cbe5f7e56272f39 100644 (file)
@@ -234,7 +234,7 @@ namespace FETools
   /**
    * This is a rather specialized function used during the construction of
    * finite element objects. It is used to build the basis of shape functions
-   * for an element given a set of polynomials and interpolation points. The
+   * for an element, given a set of polynomials and interpolation points. The
    * function is only implemented for finite elements with exactly @p dim
    * vector components. In particular, this applies to classes derived from
    * the FE_PolyTensor class.
@@ -291,7 +291,8 @@ namespace FETools
    *   that returns a scalar.
    * - That the finite element has exactly @p dim vector components.
    * - That the function $f_j$ is given by whatever the element implements
-   *   through the FiniteElement::interpolate() function.
+   *   through the FiniteElement::convert_generalized_support_point_values_to_nodal_values()
+   *   function.
    *
    * @param fe The finite element for which the operations above are to be
    *        performed.
index 03b75d83672595dff4dc19d623837633f6245f17..438fa43bdc72f47a9a1f489cf0fda23229201d57 100644 (file)
@@ -1627,8 +1627,8 @@ namespace FETools
     // We need the values of the polynomials in all generalized support points.
     // This function specifically works for the case where shape functions
     // have 'dim' vector components, so allocate that much space
-    std::vector<std::vector<double> >
-    values (dim, std::vector<double>(points.size()));
+    std::vector<Vector<double> >
+    values (points.size(), Vector<double>(dim));
 
     // In this vector, we store the
     // result of the interpolation
@@ -1645,11 +1645,12 @@ namespace FETools
         for (unsigned int k=0; k<points.size(); ++k)
           for (unsigned int d=0; d<dim; ++d)
             {
-              values[d][k] = fe.shape_value_component(i, points[k], d);
-              Assert (numbers::is_finite(values[d][k]), ExcInternalError());
+              values[k][d] = fe.shape_value_component(i, points[k], d);
+              Assert (numbers::is_finite(values[k][d]), ExcInternalError());
             }
 
-        fe.interpolate(local_dofs, values);
+        fe.convert_generalized_support_point_values_to_nodal_values(values,
+                                                                    local_dofs);
 
         // Enter the interpolated dofs into the matrix
         for (unsigned int j=0; j<n_dofs; ++j)

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.