]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix missing instantiations. 388/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 30 Dec 2014 00:32:48 +0000 (18:32 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 30 Dec 2014 00:32:48 +0000 (18:32 -0600)
include/deal.II/numerics/vector_tools.templates.h
source/numerics/vector_tools_interpolate.inst.in

index 3ecba373b02a03f6c93f1965f96377d66bb2da31..188d615ad39ff3854e88c57f78cfd16ae2075e80 100644 (file)
@@ -602,7 +602,7 @@ namespace VectorTools
     Assert(u2.size()==dof2.n_dofs(),
            ExcDimensionMismatch(u2.size(),dof2.n_dofs()));
 
-    Vector cache;
+    dealii::Vector<typename Vector::value_type> cache;
 
     // Looping over the finest common
     // mesh, this means that source and
index 1da0f724dde42cd197ad5a6687f1b368724177f9..5430ee6abceb742fcab419b28f18eb06e08a3c53 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1998 - 2013 by the deal.II authors
+// Copyright (C) 1998 - 2014 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -85,53 +85,30 @@ for (VEC : SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS; deal_II_space_dimens
   }
 
 
-for (deal_II_dimension : DIMENSIONS)
+for (VEC : SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS)
 {
   namespace VectorTools \{
 
   template
     void interpolate_to_different_mesh
     (const DoFHandler<deal_II_dimension> &,
-     const Vector<float>                 &,
+     const VEC                           &,
      const DoFHandler<deal_II_dimension> &,
-     Vector<float>                       &);
+     VEC                                 &);
 
   template
     void interpolate_to_different_mesh
     (const DoFHandler<deal_II_dimension> &,
-     const Vector<float>                 &,
+     const VEC                           &,
      const DoFHandler<deal_II_dimension> &,
      const ConstraintMatrix              &,
-     Vector<float>                       &);
+     VEC                                 &);
 
   template
     void interpolate_to_different_mesh
     (const InterGridMap<DoFHandler<deal_II_dimension> > &,
-     const Vector<float>                                &,
+     const VEC                                          &,
      const ConstraintMatrix                             &,
-     Vector<float>                                      &);
-
-  template
-    void interpolate_to_different_mesh
-    (const DoFHandler<deal_II_dimension> &,
-     const Vector<double>                &,
-     const DoFHandler<deal_II_dimension> &,
-     Vector<double>                      &);
-
-  template
-    void interpolate_to_different_mesh
-    (const DoFHandler<deal_II_dimension> &,
-     const Vector<double>                &,
-     const DoFHandler<deal_II_dimension> &,
-     const ConstraintMatrix              &,
-     Vector<double>                      &);
-
-  template
-    void interpolate_to_different_mesh
-    (const InterGridMap<DoFHandler<deal_II_dimension> > &,
-     const Vector<double>                               &,
-     const ConstraintMatrix                             &,
-     Vector<double>                                     &);
-
+     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.