]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Convert FEValuesView::Scalar to allow for generic scalar types.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 22 Feb 2015 22:14:18 +0000 (16:14 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 8 May 2015 16:19:50 +0000 (11:19 -0500)
include/deal.II/fe/fe_values.h
source/fe/fe_values.cc
source/fe/fe_values.decl.1.inst.in
source/fe/fe_values.decl.2.inst.in
source/fe/fe_values.impl.1.inst.in
source/fe/fe_values.impl.2.inst.in
source/fe/fe_values.inst.in

index 11e9e7e737f4b5fbda53500e0917e42345c7bc00..0c55d046cd577e218edb37720a384b0154dc1cf6 100644 (file)
@@ -264,11 +264,17 @@ namespace FEValuesViews
      * FEValuesBase::get_function_values function but it only works on the
      * selected scalar component.
      *
+     * The data type stored by the output vector must be what you get
+     * when you multiply the values of shape function (i.e., @p
+     * value_type) times the type used to store the values of the
+     * unknowns $U_j$ of your finite element vector $U$ (represented
+     * by the @p fe_function argument).
+     *
      * @dealiiRequiresUpdateFlags{update_values}
      */
     template <class InputVector>
     void get_function_values (const InputVector &fe_function,
-                              std::vector<value_type> &values) const;
+                              std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &values) const;
 
     /**
      * Return the gradients of the selected scalar component of the finite
@@ -280,11 +286,17 @@ namespace FEValuesViews
      * FEValuesBase::get_function_gradients function but it only works on the
      * selected scalar component.
      *
+     * The data type stored by the output vector must be what you get
+     * when you multiply the gradients of shape function (i.e., @p
+     * gradient_type) times the type used to store the values of the
+     * unknowns $U_j$ of your finite element vector $U$ (represented
+     * by the @p fe_function argument).
+     *
      * @dealiiRequiresUpdateFlags{update_gradients}
      */
     template <class InputVector>
     void get_function_gradients (const InputVector &fe_function,
-                                 std::vector<gradient_type> &gradients) const;
+                                 std::vector<typename ProductType<gradient_type,typename InputVector::value_type>::type> &gradients) const;
 
     /**
      * Return the Hessians of the selected scalar component of the finite
@@ -296,11 +308,17 @@ namespace FEValuesViews
      * FEValuesBase::get_function_hessians function but it only works on the
      * selected scalar component.
      *
+     * The data type stored by the output vector must be what you get
+     * when you multiply the Hessians of shape function (i.e., @p
+     * hessian_type) times the type used to store the values of the
+     * unknowns $U_j$ of your finite element vector $U$ (represented
+     * by the @p fe_function argument).
+     *
      * @dealiiRequiresUpdateFlags{update_hessians}
      */
     template <class InputVector>
     void get_function_hessians (const InputVector &fe_function,
-                                std::vector<hessian_type> &hessians) const;
+                                std::vector<typename ProductType<hessian_type,typename InputVector::value_type>::type> &hessians) const;
 
     /**
      * Return the Laplacians of the selected scalar component of the finite
@@ -313,11 +331,17 @@ namespace FEValuesViews
      * FEValuesBase::get_function_laplacians function but it only works on the
      * selected scalar component.
      *
+     * The data type stored by the output vector must be what you get
+     * when you multiply the Laplacians of shape function (i.e., @p
+     * value_type) times the type used to store the values of the
+     * unknowns $U_j$ of your finite element vector $U$ (represented
+     * by the @p fe_function argument).
+     *
      * @dealiiRequiresUpdateFlags{update_hessians}
      */
     template <class InputVector>
     void get_function_laplacians (const InputVector &fe_function,
-                                  std::vector<value_type> &laplacians) const;
+                                  std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &laplacians) const;
 
   private:
     /**
index 858f50a69e3bd27ee54b259b35a5c77cd112d6da..04536aff69a79dd887bd280e820797c4816ebe51 100644 (file)
@@ -445,22 +445,16 @@ namespace FEValuesViews
 
   namespace internal
   {
-    // put the evaluation part of the get_function_xxx from a local vector
-    // into separate functions. this reduces the size of the compilation unit
-    // by a factor more than 2 without affecting the performance at all.
-
-    // remark: up to revision 27774, dof_values used to be extracted as
-    // VectorType::value_type and not simply double. this did not make a lot
-    // of sense since they were later extracted and converted to double
-    // consistently throughout the code since revision 17903 at least.
+    // Given values of degrees of freedom, evaluate the
+    // values/gradients/... at quadrature points
 
     // ------------------------- scalar functions --------------------------
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_values (const ::dealii::Vector<double> &dof_values,
+    do_function_values (const ::dealii::Vector<Number> &dof_values,
                         const Table<2,double>          &shape_values,
                         const std::vector<typename Scalar<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                        std::vector<double>            &values)
+                        std::vector<typename ProductType<Number,double>::type>            &values)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -488,12 +482,12 @@ namespace FEValuesViews
 
     // same code for gradient and Hessian, template argument 'order' to give
     // the order of the derivative (= rank of gradient/Hessian tensor)
-    template <int order, int dim, int spacedim>
+    template <int order, int dim, int spacedim, typename Number>
     void
-    do_function_derivatives (const ::dealii::Vector<double> &dof_values,
+    do_function_derivatives (const ::dealii::Vector<Number> &dof_values,
                              const std::vector<std::vector<dealii::Tensor<order,spacedim> > > &shape_derivatives,
                              const std::vector<typename Scalar<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                             std::vector<dealii::Tensor<order,spacedim> > &derivatives)
+                             std::vector<typename ProductType<Number,dealii::Tensor<order,spacedim> >::type> &derivatives)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -514,18 +508,19 @@ namespace FEValuesViews
             const dealii::Tensor<order,spacedim> *shape_derivative_ptr =
               &shape_derivatives[shape_function_data[shape_function].row_index][0];
             for (unsigned int q_point=0; q_point<n_quadrature_points; ++q_point)
-              derivatives[q_point] += value **shape_derivative_ptr++;
+              derivatives[q_point] += value *
+                                     typename ProductType<Number,dealii::Tensor<order,spacedim> >::type(*shape_derivative_ptr++);
           }
     }
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_laplacians (const ::dealii::Vector<double> &dof_values,
+    do_function_laplacians (const ::dealii::Vector<Number> &dof_values,
                             const std::vector<std::vector<dealii::Tensor<2,spacedim> > > &shape_hessians,
                             const std::vector<typename Scalar<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                            std::vector<double>           &laplacians)
+                            std::vector<typename ProductType<Number,double>::type>           &laplacians)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -553,11 +548,11 @@ namespace FEValuesViews
 
     // ----------------------------- vector part ---------------------------
 
-    template <int dim, int spacedim>
-    void do_function_values (const ::dealii::Vector<double> &dof_values,
+    template <int dim, int spacedim, typename Number>
+    void do_function_values (const ::dealii::Vector<Number> &dof_values,
                              const Table<2,double>          &shape_values,
                              const std::vector<typename Vector<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                             std::vector<dealii::Tensor<1,spacedim> > &values)
+                             std::vector<typename ProductType<Number,dealii::Tensor<1,spacedim> >::type> &values)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -601,12 +596,12 @@ namespace FEValuesViews
 
 
 
-    template <int order, int dim, int spacedim>
+    template <int order, int dim, int spacedim, typename Number>
     void
-    do_function_derivatives (const ::dealii::Vector<double> &dof_values,
+    do_function_derivatives (const ::dealii::Vector<Number> &dof_values,
                              const std::vector<std::vector<dealii::Tensor<order,spacedim> > > &shape_derivatives,
                              const std::vector<typename Vector<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                             std::vector<dealii::Tensor<order+1,spacedim> > &derivatives)
+                             std::vector<typename ProductType<Number,dealii::Tensor<order+1,spacedim> >::type> &derivatives)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -653,12 +648,12 @@ namespace FEValuesViews
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_symmetric_gradients (const ::dealii::Vector<double> &dof_values,
+    do_function_symmetric_gradients (const ::dealii::Vector<Number> &dof_values,
                                      const std::vector<std::vector<dealii::Tensor<1,spacedim> > > &shape_gradients,
                                      const std::vector<typename Vector<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                                     std::vector<dealii::SymmetricTensor<2,spacedim> > &symmetric_gradients)
+                                     std::vector<typename ProductType<Number,dealii::SymmetricTensor<2,spacedim> >::type> &symmetric_gradients)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -706,12 +701,12 @@ namespace FEValuesViews
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_divergences (const ::dealii::Vector<double> &dof_values,
+    do_function_divergences (const ::dealii::Vector<Number> &dof_values,
                              const std::vector<std::vector<dealii::Tensor<1,spacedim> > > &shape_gradients,
                              const std::vector<typename Vector<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                             std::vector<double> &divergences)
+                             std::vector<typename ProductType<Number,double>::type> &divergences)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -756,12 +751,12 @@ namespace FEValuesViews
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_curls (const ::dealii::Vector<double> &dof_values,
+    do_function_curls (const ::dealii::Vector<Number> &dof_values,
                        const std::vector<std::vector<dealii::Tensor<1,spacedim> > > &shape_gradients,
                        const std::vector<typename Vector<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                       std::vector<typename dealii::internal::CurlType<spacedim>::type> &curls)
+                       std::vector<typename ProductType<Number,typename dealii::internal::CurlType<spacedim>::type>::type> &curls)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -948,12 +943,12 @@ namespace FEValuesViews
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_laplacians (const ::dealii::Vector<double> &dof_values,
+    do_function_laplacians (const ::dealii::Vector<Number> &dof_values,
                             const std::vector<std::vector<dealii::Tensor<2,spacedim> > > &shape_hessians,
                             const std::vector<typename Vector<dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                            std::vector<dealii::Tensor<1,spacedim> > &laplacians)
+                            std::vector<typename ProductType<Number,dealii::Tensor<1,spacedim> >::type> &laplacians)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -1002,12 +997,12 @@ namespace FEValuesViews
 
     // ---------------------- symmetric tensor part ------------------------
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_values (const ::dealii::Vector<double> &dof_values,
+    do_function_values (const ::dealii::Vector<Number> &dof_values,
                         const dealii::Table<2,double>          &shape_values,
                         const std::vector<typename SymmetricTensor<2,dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                        std::vector<dealii::SymmetricTensor<2,spacedim> > &values)
+                        std::vector<typename ProductType<Number,dealii::SymmetricTensor<2,spacedim> >::type> &values)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -1056,12 +1051,12 @@ namespace FEValuesViews
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_divergences (const ::dealii::Vector<double> &dof_values,
+    do_function_divergences (const ::dealii::Vector<Number> &dof_values,
                              const std::vector<std::vector<dealii::Tensor<1,spacedim> > > &shape_gradients,
                              const std::vector<typename SymmetricTensor<2,dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                             std::vector<dealii::Tensor<1,spacedim> > &divergences)
+                             std::vector<typename ProductType<Number,dealii::Tensor<1,spacedim> >::type> &divergences)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -1144,12 +1139,12 @@ namespace FEValuesViews
 
     // ---------------------- non-symmetric tensor part ------------------------
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_values (const ::dealii::Vector<double> &dof_values,
+    do_function_values (const ::dealii::Vector<Number> &dof_values,
                         const dealii::Table<2,double>          &shape_values,
                         const std::vector<typename Tensor<2,dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                        std::vector<dealii::Tensor<2,spacedim> > &values)
+                        std::vector<typename ProductType<Number,dealii::Tensor<2,spacedim> >::type> &values)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -1200,12 +1195,12 @@ namespace FEValuesViews
 
 
 
-    template <int dim, int spacedim>
+    template <int dim, int spacedim, typename Number>
     void
-    do_function_divergences (const ::dealii::Vector<double> &dof_values,
+    do_function_divergences (const ::dealii::Vector<Number> &dof_values,
                              const std::vector<std::vector<dealii::Tensor<1,spacedim> > > &shape_gradients,
                              const std::vector<typename Tensor<2,dim,spacedim>::ShapeFunctionData> &shape_function_data,
-                             std::vector<dealii::Tensor<1,spacedim> > &divergences)
+                             std::vector<typename ProductType<Number,dealii::Tensor<1,spacedim> >::type> &divergences)
     {
       const unsigned int dofs_per_cell = dof_values.size();
       const unsigned int n_quadrature_points = dofs_per_cell > 0 ?
@@ -1267,7 +1262,7 @@ namespace FEValuesViews
   void
   Scalar<dim,spacedim>::
   get_function_values (const InputVector &fe_function,
-                       std::vector<value_type> &values) const
+                       std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &values) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
     Assert (fe_values.update_flags & update_values,
@@ -1278,7 +1273,7 @@ namespace FEValuesViews
                      fe_values.present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell and call internal worker function
-    dealii::Vector<double> dof_values(fe_values.dofs_per_cell);
+    dealii::Vector<typename InputVector::value_type> dof_values(fe_values.dofs_per_cell);
     fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_values<dim,spacedim>
     (dof_values, fe_values.shape_values, shape_function_data, values);
@@ -1291,7 +1286,7 @@ namespace FEValuesViews
   void
   Scalar<dim,spacedim>::
   get_function_gradients (const InputVector &fe_function,
-                          std::vector<gradient_type> &gradients) const
+                          std::vector<typename ProductType<gradient_type,typename InputVector::value_type>::type> &gradients) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
     Assert (fe_values.update_flags & update_gradients,
@@ -1302,7 +1297,7 @@ namespace FEValuesViews
                      fe_values.present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<double> dof_values (fe_values.dofs_per_cell);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
     fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<1,dim,spacedim>
     (dof_values, fe_values.shape_gradients, shape_function_data, gradients);
@@ -1315,7 +1310,7 @@ namespace FEValuesViews
   void
   Scalar<dim,spacedim>::
   get_function_hessians (const InputVector &fe_function,
-                         std::vector<hessian_type> &hessians) const
+                         std::vector<typename ProductType<hessian_type,typename InputVector::value_type>::type> &hessians) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
     Assert (fe_values.update_flags & update_hessians,
@@ -1326,7 +1321,7 @@ namespace FEValuesViews
                      fe_values.present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<double> dof_values (fe_values.dofs_per_cell);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
     fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<2,dim,spacedim>
     (dof_values, fe_values.shape_hessians, shape_function_data, hessians);
@@ -1339,7 +1334,7 @@ namespace FEValuesViews
   void
   Scalar<dim,spacedim>::
   get_function_laplacians (const InputVector &fe_function,
-                           std::vector<value_type> &laplacians) const
+                           std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &laplacians) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
     Assert (fe_values.update_flags & update_hessians,
@@ -1350,7 +1345,7 @@ namespace FEValuesViews
                      fe_values.present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<double> dof_values (fe_values.dofs_per_cell);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
     fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_laplacians<dim,spacedim>
     (dof_values, fe_values.shape_hessians, shape_function_data, laplacians);
@@ -1769,7 +1764,7 @@ public:
   virtual
   void
   get_interpolated_dof_values (const IndexSet &in,
-                               Vector<double> &out) const = 0;
+                               Vector<IndexSet::value_type> &out) const = 0;
 };
 
 
@@ -1838,7 +1833,7 @@ public:
   virtual
   void
   get_interpolated_dof_values (const IndexSet &in,
-                               Vector<double> &out) const;
+                               Vector<IndexSet::value_type> &out) const;
 
 private:
   /**
@@ -1941,7 +1936,7 @@ public:
   virtual
   void
   get_interpolated_dof_values (const IndexSet &in,
-                               Vector<double> &out) const;
+                               Vector<IndexSet::value_type> &out) const;
 
 private:
   /**
@@ -2005,7 +2000,7 @@ template <typename CI>
 void
 FEValuesBase<dim,spacedim>::CellIterator<CI>::
 get_interpolated_dof_values (const IndexSet &in,
-                             Vector<double> &out) const
+                             Vector<IndexSet::value_type> &out) const
 {
   Assert (cell->has_children() == false, ExcNotImplemented());
 
@@ -2065,7 +2060,7 @@ template <int dim, int spacedim>
 void
 FEValuesBase<dim,spacedim>::TriaCellIterator::
 get_interpolated_dof_values (const IndexSet &,
-                             Vector<double> &) const
+                             Vector<IndexSet::value_type> &) const
 {
   Assert (false, ExcMessage (message_string));
 }
index b27b8380166e31bc57387aae8d9f08c305d86200..30129eb20a3c104696802382fe11d499e471dbd6 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2010 - 2014 by the deal.II authors
+// Copyright (C) 2010 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -27,5 +27,5 @@ for (VEC : SERIAL_VECTORS)
     virtual
     void
     get_interpolated_dof_values (const VEC &in,
-                                Vector<double> &out) const = 0;
+                                Vector<VEC::value_type> &out) const = 0;
   }
index d0ba5103d478c1d0ae2c09bd707116c3699f3c1b..9a91f9830fb22d6249ee3182bbbba6d73b3fb1ef 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2010 - 2014 by the deal.II authors
+// Copyright (C) 2010 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -27,5 +27,5 @@ for (VEC : SERIAL_VECTORS)
     virtual
     void
     get_interpolated_dof_values (const VEC      &in,
-                                Vector<double> &out) const;
+                                Vector<VEC::value_type> &out) const;
   }
index abeca627c63da0334d45916085889a7a2318ee6a..7dbdc341cf4ae66a6de5dd9a809366c09f9b0c80 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2010 - 2014 by the deal.II authors
+// Copyright (C) 2010 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -21,7 +21,7 @@ for (VEC : SERIAL_VECTORS)
     void
     FEValuesBase<dim,spacedim>::CellIterator<CI>::
     get_interpolated_dof_values (const VEC      &in,
-                                Vector<double> &out) const
+                                Vector<VEC::value_type> &out) const
     \{
       cell->get_interpolated_dof_values (in, out);
     \}
index 363bff5a850a4f6c54bbdfff816875dfa205a9a0..5aa37dbc4d9dde0c014348be7dcf4116ff391707 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2010 - 2014 by the deal.II authors
+// Copyright (C) 2010 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -20,7 +20,7 @@ for (VEC : SERIAL_VECTORS)
     void
     FEValuesBase<dim,spacedim>::TriaCellIterator::
     get_interpolated_dof_values (const VEC &,
-                                Vector<double> &) const
+                                Vector<VEC::value_type> &) const
     \{
       Assert (false, ExcMessage (message_string));
     \}
index fb0624352163155ee6941d068954f19f12e3ee29..61716fa11f03a7e66c0ac28b4e9ad73cf099d456 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1998 - 2014 by the deal.II authors
+// Copyright (C) 1998 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -70,19 +70,19 @@ for (VEC : SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimensi
     template
       void FEValuesViews::Scalar<deal_II_dimension, deal_II_space_dimension>
       ::get_function_values<dealii::VEC>
-      (const dealii::VEC&, std::vector<value_type>&) const;
+      (const dealii::VEC&, std::vector<typename ProductType<typename dealii::VEC::value_type,value_type>::type>&) const;
     template
       void FEValuesViews::Scalar<deal_II_dimension, deal_II_space_dimension>
       ::get_function_gradients<dealii::VEC>
-      (const dealii::VEC&, std::vector<dealii::Tensor<1,deal_II_space_dimension> >&) const;
+      (const dealii::VEC&, std::vector<typename ProductType<typename dealii::VEC::value_type,dealii::Tensor<1,deal_II_space_dimension>>::type>&) const;
     template
       void FEValuesViews::Scalar<deal_II_dimension, deal_II_space_dimension>
       ::get_function_hessians<dealii::VEC>
-      (const dealii::VEC&, std::vector<dealii::Tensor<2,deal_II_space_dimension> >&) const;
+      (const dealii::VEC&, std::vector<typename ProductType<typename dealii::VEC::value_type,dealii::Tensor<2,deal_II_space_dimension>>::type>&) const;
     template
       void FEValuesViews::Scalar<deal_II_dimension, deal_II_space_dimension>
       ::get_function_laplacians<dealii::VEC>
-      (const dealii::VEC&, std::vector<double>&) const;
+      (const dealii::VEC&, std::vector<typename ProductType<typename dealii::VEC::value_type,double>::type> &) const;
 
     template
       void FEValuesViews::Vector<deal_II_dimension, deal_II_space_dimension>
@@ -261,16 +261,16 @@ for (deal_II_dimension : DIMENSIONS)
 #if (deal_II_dimension == DIM_A) || (deal_II_dimension == DIM_B)
     template
        void FEValuesViews::Scalar<deal_II_dimension>::get_function_values<dealii::IndexSet>
-       (const dealii::IndexSet&, std::vector<double>&) const;
+       (const dealii::IndexSet&, std::vector<typename ProductType<IndexSet::value_type,double>::type> &) const;
     template
        void FEValuesViews::Scalar<deal_II_dimension>::get_function_gradients<dealii::IndexSet>
-       (const dealii::IndexSet&, std::vector<dealii::Tensor<1,deal_II_dimension> >&) const;
+       (const dealii::IndexSet&, std::vector<typename ProductType<IndexSet::value_type,dealii::Tensor<1,deal_II_dimension> >::type>&) const;
     template
        void FEValuesViews::Scalar<deal_II_dimension>::get_function_hessians<dealii::IndexSet>
-       (const dealii::IndexSet&, std::vector<dealii::Tensor<2,deal_II_dimension> >&) const;
+       (const dealii::IndexSet&, std::vector<typename ProductType<IndexSet::value_type,dealii::Tensor<2,deal_II_dimension> >::type>&) const;
     template
        void FEValuesViews::Scalar<deal_II_dimension>::get_function_laplacians<dealii::IndexSet>
-       (const dealii::IndexSet&, std::vector<double>&) const;
+       (const dealii::IndexSet&, std::vector<typename ProductType<IndexSet::value_type,double>::type>&) const;
 
     template
        void FEValuesViews::Vector<deal_II_dimension>::get_function_values<dealii::IndexSet>

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.