]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Work around the Intel ICC bug that leads to forgotten template instantiations.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 27 Nov 2013 15:17:29 +0000 (15:17 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 27 Nov 2013 15:17:29 +0000 (15:17 +0000)
git-svn-id: https://svn.dealii.org/trunk@31811 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/numerics/vector_tools.h
deal.II/include/deal.II/numerics/vector_tools.templates.h
deal.II/source/numerics/vector_tools_interpolate.inst.in

index 89ec0b712c427a3ed3d007049be3934ba7ae2064..26b054bb08d3742b5d22523459d55326b83f4afa 100644 (file)
@@ -40,6 +40,19 @@ inconvenience this causes.
 </p>
 
 <ol>
+  <li>
+  Changed: The kinds of template arguments for the VectorTools::interpolate
+  function taking a Mapping as first argument has changed. This was done to
+  work around a bug in the Intel ICC compiler which led to linker errors. Since
+  the actual function argument list remains unchanged, the only way you will
+  notice this change is if you <i>explicitly</i> specify template arguments.
+  The only place one would typically do that is if you take the address of
+  a template function. Since this is not a common operation, the impact of this
+  change is probably limited.
+  <br>
+  (Wolfgang Bangerth, 2013/11/27)
+  </li>
+
   <li>
   Changed: The ghost handling of the parallel::distributed::Vector class has
   been reworked: The vector now carries a global state that stores whether
index daf77a045dfadcb705757542527a0ddab47f6924..28c5fa3f10fe69cd4e0702a9ecade5fdb94784d2 100644 (file)
@@ -451,11 +451,11 @@ namespace VectorTools
    * replaced by a hp::MappingCollection in
    * case of a hp::DoFHandler.
    */
-  template <class VECTOR, class DH>
-  void interpolate (const Mapping<DH::dimension,DH::space_dimension>    &mapping,
-                    const DH              &dof,
-                    const Function<DH::space_dimension>   &function,
-                    VECTOR                &vec);
+  template <class VECTOR, int dim, int spacedim, template <int,int> class DH>
+  void interpolate (const Mapping<dim,spacedim>    &mapping,
+                    const DH<dim,spacedim>         &dof,
+                    const Function<spacedim>       &function,
+                    VECTOR                         &vec);
 
   /**
    * Calls the @p interpolate()
index a3748c86b89a4210f2e1643058f12fd076e5307d..77cc717739a52e512e661f0f3a1a74888d1f9569 100644 (file)
@@ -74,24 +74,22 @@ DEAL_II_NAMESPACE_OPEN
 namespace VectorTools
 {
 
-  template <class VECTOR, class DH>
-  void interpolate (const Mapping<DH::dimension,DH::space_dimension>    &mapping,
-                    const DH              &dof,
-                    const Function<DH::space_dimension>   &function,
-                    VECTOR                &vec)
+  template <class VECTOR, int dim, int spacedim, template <int,int> class DH>
+  void interpolate (const Mapping<dim,spacedim>    &mapping,
+                    const DH<dim,spacedim>         &dof,
+                    const Function<spacedim>       &function,
+                    VECTOR                         &vec)
   {
-    const unsigned int dim=DH::dimension;
-
     Assert (dof.get_fe().n_components() == function.n_components,
             ExcDimensionMismatch(dof.get_fe().n_components(),
                                  function.n_components));
 
-    const hp::FECollection<DH::dimension,DH::space_dimension> fe (dof.get_fe());
+    const hp::FECollection<dim,spacedim> fe (dof.get_fe());
     const unsigned int          n_components = fe.n_components();
     const bool                  fe_is_system = (n_components != 1);
 
-    typename DH::active_cell_iterator cell = dof.begin_active(),
-                                      endc = dof.end();
+    typename DH<dim,spacedim>::active_cell_iterator cell = dof.begin_active(),
+                                                   endc = dof.end();
 
     // For FESystems many of the
     // unit_support_points will appear
@@ -180,7 +178,7 @@ namespace VectorTools
     const unsigned int max_rep_points = *std::max_element (n_rep_points.begin(),
                                                            n_rep_points.end());
     std::vector<types::global_dof_index> dofs_on_cell (fe.max_dofs_per_cell());
-    std::vector<Point<DH::space_dimension> >  rep_points (max_rep_points);
+    std::vector<Point<spacedim> >  rep_points (max_rep_points);
 
     // get space for the values of the
     // function at the rep support points.
@@ -199,9 +197,9 @@ namespace VectorTools
 
     // Transformed support points are computed by
     // FEValues
-    hp::MappingCollection<dim,DH::space_dimension> mapping_collection (mapping);
+    hp::MappingCollection<dim,spacedim> mapping_collection (mapping);
 
-    hp::FEValues<dim, DH::space_dimension> fe_values (mapping_collection,
+    hp::FEValues<dim,spacedim> fe_values (mapping_collection,
                                                       fe, support_quadrature, update_quadrature_points);
 
     for (; cell!=endc; ++cell)
@@ -213,7 +211,7 @@ namespace VectorTools
           // get location of finite element
           // support_points
           fe_values.reinit(cell);
-          const std::vector<Point<DH::space_dimension> > &support_points =
+          const std::vector<Point<spacedim> > &support_points =
             fe_values.get_present_fe_values().get_quadrature_points();
 
           // pick out the representative
@@ -5218,11 +5216,11 @@ namespace VectorTools
         case H1_norm:
           exponent = 2.;
           break;
-         
+
         case L1_norm:
           exponent = 1.;
           break;
-         
+
         default:
           break;
         }
@@ -5247,7 +5245,7 @@ namespace VectorTools
           if (spacedim == dim+1)
            update_flags |= UpdateFlags (update_normal_vectors);
           // no break!
-         
+
         default:
           update_flags |= UpdateFlags (update_values);
           break;
index e44dd85c4ec5dcae13d20c2d81e6f92f3d69650d..bbc0d6b237bbe2ce3ed1ffdfbe33c8ee07a2372e 100644 (file)
@@ -22,7 +22,7 @@ for (VEC : SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS; deal_II_space_dimens
 
       template
         void interpolate
-        (const Mapping<DoFHandler<deal_II_dimension,deal_II_space_dimension>::dimension,DoFHandler<deal_II_dimension,deal_II_space_dimension>::space_dimension>&,
+        (const Mapping<deal_II_dimension,deal_II_space_dimension> &,
          const DoFHandler<deal_II_dimension,deal_II_space_dimension>&,
          const Function<deal_II_space_dimension>&,
          VEC&);

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.